Peano
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
AddVolumeAndFaceSolution.py
Go to the documentation of this file.
1# This file is part of the ExaHyPE2 project. For conditions of distribution and
2# use, please see the copyright notice at www.peano-framework.org
3
4from .AbstractRungeKuttaDGActionSet import AbstractRungeKuttaDGActionSet
5
6from exahype2.solvers.ButcherTableau import ButcherTableau
7
8import peano4
9import jinja2
10
12 """
13
14 The linear combination of the Runge Kutta trials has to be projected onto
15 the faces, so we can then solve the Riemann problems. So the projection
16 happens in one grid sweep, the corresponding Riemann solve in the next one.
17
18
19 """
20
21 TemplateProject = """
22 {% for PREDICATE_NO in range(0,PREDICATES|length) %}
23 if ({{PREDICATES[PREDICATE_NO]}}) {
24 const double timeStampOldSolution = fineGridCell{{SOLVER_NAME}}CellLabel.getTimeStamp();
25
26 // Set the variable
27 // double timeStepSize
28 {{COMPUTE_TIME_STEP_SIZE}}
29
30 const double timeStamp = timeStampOldSolution + {{BUTCHER_TABLEAU_RELATIVE_TIME_STEP_SIZES[PREDICATE_NO]}} * timeStepSize;
31
32 assertion3( tarch::la::greaterEquals( timeStepSize, 0.0 ), timeStamp, timeStepSize, timeStampOldSolution );
33 assertion3( tarch::la::greaterEquals( timeStamp, 0.0 ), timeStamp, timeStepSize, timeStampOldSolution );
34
35 #if Dimensions==2
36 double* dQdt = fineGridCell{{UNKNOWN_IDENTIFIER}}RhsEstimates.value + {{PREDICATE_NO}} * {{NUMBER_OF_DOFS_PER_CELL_2D}} * {{NUMBER_OF_UNKNOWNS}};
37 #elif Dimensions==3
38 double* dQdt = fineGridCell{{UNKNOWN_IDENTIFIER}}RhsEstimates.value + {{PREDICATE_NO}} * {{NUMBER_OF_DOFS_PER_CELL_3D}} * {{NUMBER_OF_UNKNOWNS}};
39 #endif
40
41 ::exahype2::dg::{{KERNEL_NAMESPACE}}::{{ADD_SOLVER_CONTRIBUTIONS_CALL}}
42
43 ::exahype2::CellData<double, double> cellData(
44 nullptr,
45 marker.x(),
46 marker.h(),
47 timeStamp,
48 timeStepSize,
49 dQdt
50 );
51
52 ::exahype2::dg::{{KERNEL_NAMESPACE}}::{{MULTIPLY_WITH_INVERTED_MASS_MATRIX_CALL}}
53
54 const double* oldQ = fineGridCell{{UNKNOWN_IDENTIFIER}}.value;
55
56 {{POSTPROCESS_UPDATED_CELL_AFTER_RUNGE_KUTTA_STEP}}
57 }
58 {% endfor %}
59"""
60
61
62 def __init__(self,solver):
63 """
64
65 guard_project: String (C++ code)
66 Predicate which controls if the solution is actually projected
67
68 guard_safe_old_time_step: String (C++ code)
69 Predicate which controls if the projection should be copied into
70 the old solution and the time step should also be moved over
71
72 """
73 super(AddVolumeAndFaceSolution,self).__init__(solver)
76
77
78 @property
79 def guards(self):
80 if self._guards==[]:
81 raise Exception( "Guards are not initialised" )
82 return self._guards
83
84
85 @guards.setter
86 def guards(self,new_guards):
87 if new_guards!=[] and len(new_guards)!=self._solver.number_of_Runge_Kutta_steps():
88 raise Exception( "Expect one guard per Runge Kutta step. Have {} steps but got guards {}".format(solver.number_of_Runge_Kutta_steps(),guards) )
89 self._guards = new_guards
90
91
92 def get_body_of_operation(self,operation_name):
93 result = ""
94 if operation_name==peano4.solversteps.ActionSet.OPERATION_TOUCH_CELL_FIRST_TIME:
95 d = {}
96 self._solver._init_dictionary_with_default_parameters(d)
97 self._solver.add_entries_to_text_replacement_dictionary(d)
98 d[ "PREDICATES" ] = self._guards
99 d[ "BUTCHER_TABLEAU_WEIGHTS" ] = self._butcher_tableau.weight_matrix()
100 d[ "BUTCHER_TABLEAU_RELATIVE_TIME_STEP_SIZES" ] = self._butcher_tableau.time_step_sizes()
101 result = jinja2.Template(self.TemplateProject).render(**d)
102 pass
103 return result
104
105
107 return __name__.replace(".py", "").replace(".", "_")
The linear combination of the Runge Kutta trials has to be projected onto the faces,...
__init__(self, solver)
guard_project: String (C++ code) Predicate which controls if the solution is actually projected
get_action_set_name(self)
You should replicate this function in each subclass, so you get meaningful action set names (otherwis...
get_body_of_operation(self, operation_name)
Return actual C++ code snippets to be inserted into C++ code.