From ffc071d2bc9f80330e0357a31b9dfe54325dcc75 Mon Sep 17 00:00:00 2001 From: =?UTF-8?q?Ren=C3=A9=20He=C3=9F?= <rene.hess@iwr.uni-heidelberg.de> Date: Thu, 18 Aug 2016 14:16:20 +0200 Subject: [PATCH] Allow multiple forms: The driver solves the PDE given by the first form in forms. This is a first step towards instationary problems. --- python/dune/perftool/compile.py | 62 ++++++++++++-------- python/dune/perftool/pdelab/driver.py | 8 ++- python/dune/perftool/pdelab/localoperator.py | 17 +++++- test/CMakeLists.txt | 3 +- test/heatequation/CMakeLists.txt | 14 +++++ test/heatequation/heatequation.mini | 11 ++++ test/heatequation/heatequation.ufl | 11 ++++ 7 files changed, 95 insertions(+), 31 deletions(-) create mode 100644 test/heatequation/CMakeLists.txt create mode 100644 test/heatequation/heatequation.mini create mode 100644 test/heatequation/heatequation.ufl diff --git a/python/dune/perftool/compile.py b/python/dune/perftool/compile.py index 56ad0096..2284e14a 100644 --- a/python/dune/perftool/compile.py +++ b/python/dune/perftool/compile.py @@ -4,6 +4,7 @@ The methods to run the parts of the form compiler Should also contain the entrypoint methods. """ from __future__ import absolute_import +from os.path import splitext # Disable loopy caching before we do anything else! import loopy @@ -32,24 +33,25 @@ def read_ufl(uflfile): # Extract the form from the given data data = interpret_ufl_namespace(namespace) - form = data.forms[0] - formdata = compute_form_data(form) - form = formdata.preprocessed_form + formdatas = [] + forms = data.forms + for index, form in enumerate(forms): + formdatas.append(compute_form_data(form)) + forms[index] = formdatas[index].preprocessed_form - # We do not expect more than one form - assert len(data.forms) == 1 + # We expect at least one form + assert len(data.forms) >= 1 # apply some transformations unconditionally! from dune.perftool.ufl.transformations import transform_form from dune.perftool.ufl.transformations.indexpushdown import pushdown_indexed from dune.perftool.ufl.transformations.reindexing import reindexing + for i in range(len(forms)): + forms[i] = transform_form(forms[i], pushdown_indexed) + forms[i] = transform_form(forms[i], reindexing) + formdatas[i].preprocessed_form = forms[i] - form = transform_form(form, pushdown_indexed) - form = transform_form(form, reindexing) - - formdata.preprocessed_form = form - - return formdata, data + return formdatas, data def generate_driver(form, filename): @@ -79,22 +81,30 @@ def generate_driver(form, filename): # This function is the entrypoint of the ufl2pdelab executable def compile_form(): from dune.perftool.options import get_option - formdata, data = read_ufl(get_option("uflfile")) + formdatas, data = read_ufl(get_option("uflfile")) from dune.perftool.generation import cache_context, global_context - with global_context(formdata=formdata, data=data): + with global_context(data=data, formdata=formdatas[0]): + # Generate driver file if get_option("driver_file"): with cache_context('driver', delete=True): - generate_driver(formdata.preprocessed_form, get_option("driver_file")) - - if get_option("operator_file"): - from dune.perftool.pdelab.localoperator import generate_localoperator_kernels - kernels = generate_localoperator_kernels(formdata, data) - - # TODO insert sophisticated analysis/feedback loops here - if get_option("interactive"): - from dune.perftool.interactive import start_interactive_session - start_interactive_session(kernels) - - from dune.perftool.pdelab.localoperator import generate_localoperator_file - generate_localoperator_file(kernels) + generate_driver(formdatas[0].preprocessed_form, get_option("driver_file")) + + # Generate local operator files + for formdata in formdatas: + with global_context(formdata=formdata): + if get_option("operator_file"): + from dune.perftool.pdelab.localoperator import generate_localoperator_kernels + kernels = generate_localoperator_kernels(formdata, data) + + # TODO insert sophisticated analysis/feedback loops here + if get_option("interactive"): + from dune.perftool.interactive import start_interactive_session + start_interactive_session(kernels) + + if get_option("operator_file"): + from dune.perftool.pdelab.localoperator import name_localoperator_file + filename = name_localoperator_file(formdata, data) + + from dune.perftool.pdelab.localoperator import generate_localoperator_file + generate_localoperator_file(kernels, filename) diff --git a/python/dune/perftool/pdelab/driver.py b/python/dune/perftool/pdelab/driver.py index ff87bd11..d5840e6b 100644 --- a/python/dune/perftool/pdelab/driver.py +++ b/python/dune/perftool/pdelab/driver.py @@ -602,8 +602,12 @@ def typedef_localoperator(name): # return "typedef LocalOperator<{}> {};".format(params, name) ugfs = type_gfs(_form.coefficients()[0].ufl_element()) vgfs = type_gfs(_form.arguments()[0].ufl_element()) - from dune.perftool.options import get_option - include_file(get_option('operator_file'), filetag="driver") + from dune.perftool.generation import get_global_context_value + data = get_global_context_value("data") + formdata = get_global_context_value("formdata") + from dune.perftool.pdelab.localoperator import name_localoperator_file + filename = name_localoperator_file(formdata, data) + include_file(filename, filetag="driver") return "using {} = LocalOperator<{}, {}>;".format(name, ugfs, vgfs) diff --git a/python/dune/perftool/pdelab/localoperator.py b/python/dune/perftool/pdelab/localoperator.py index b539e14f..a072bc9d 100644 --- a/python/dune/perftool/pdelab/localoperator.py +++ b/python/dune/perftool/pdelab/localoperator.py @@ -1,4 +1,5 @@ from __future__ import absolute_import +from os.path import splitext from dune.perftool.options import get_option from dune.perftool.generation import (base_class, @@ -18,6 +19,18 @@ from dune.perftool.cgen.clazz import (AccessModifier, from dune.perftool import Restriction +def name_localoperator_file(formdata, data): + # Check wether the formdata has a name in UFL + try: + suffix = '_' + data.object_names[id(formdata.original_form)] + except: + suffix = '' + from dune.perftool.options import get_option + basename, extension = splitext(get_option("operator_file")) + filename = basename + suffix + extension + return filename + + @template_parameter("operator") def lop_template_ansatz_gfs(): return "GFSU" @@ -406,7 +419,7 @@ def generate_localoperator_kernels(formdata, data): return operator_kernels -def generate_localoperator_file(kernels): +def generate_localoperator_file(kernels, filename): operator_methods = [] # Make generables from the given kernels @@ -422,4 +435,4 @@ def generate_localoperator_file(kernels): param = cgen_class_from_cache("parameterclass") # TODO take the name of this thing from the UFL file lop = cgen_class_from_cache("operator", members=operator_methods) - generate_file(get_option("operator_file"), "operatorfile", [param, lop]) + generate_file(filename, "operatorfile", [param, lop]) diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt index 76611b31..d41b9394 100644 --- a/test/CMakeLists.txt +++ b/test/CMakeLists.txt @@ -1,4 +1,5 @@ +add_subdirectory(heatequation) add_subdirectory(laplace) +add_subdirectory(nonlinear) add_subdirectory(poisson) add_subdirectory(stokes) -add_subdirectory(nonlinear) diff --git a/test/heatequation/CMakeLists.txt b/test/heatequation/CMakeLists.txt new file mode 100644 index 00000000..f10d2d13 --- /dev/null +++ b/test/heatequation/CMakeLists.txt @@ -0,0 +1,14 @@ +add_generated_executable(UFLFILE heatequation.ufl + TARGET heatequation + FORM_COMPILER_ARGS + ) + + +dune_add_system_test(TARGET heatequation + INIFILE heatequation.mini + ) + +# # Add the reference code with the PDELab localoperator that produced +# # the reference vtk file +# add_executable(heatequation_ref heatequation_main.cc) +# set_target_properties(heatequation_ref PROPERTIES EXCLUDE_FROM_ALL 1) diff --git a/test/heatequation/heatequation.mini b/test/heatequation/heatequation.mini new file mode 100644 index 00000000..00a89305 --- /dev/null +++ b/test/heatequation/heatequation.mini @@ -0,0 +1,11 @@ +__name = heatequation + +lowerleft = 0.0 0.0 +upperright = 1.0 1.0 +elements = 32 32 +elementType = simplical + +[wrapper.vtkcompare] +name = {__name} +reference = heatequation_ref +extension = vtu diff --git a/test/heatequation/heatequation.ufl b/test/heatequation/heatequation.ufl new file mode 100644 index 00000000..922c3276 --- /dev/null +++ b/test/heatequation/heatequation.ufl @@ -0,0 +1,11 @@ +f = Expression("Dune::FieldVector<double,2> c(0.5); c-= x; return 4.*(1.-c.two_norm2())*std::exp(-1.*c.two_norm2());") +g = Expression("Dune::FieldVector<double,2> c(0.5); c-= x; return std::exp(-1.*c.two_norm2());") + +V = FiniteElement("CG", "triangle", 1, dirichlet_expression=g) +u = TrialFunction(V) +v = TestFunction(V) + +mass = (u*v)*dx +poisson = (inner(grad(u), grad(v)) - f*v)*dx + +forms = [poisson,mass] -- GitLab