Peano
Loading...
Searching...
No Matches
FiniteVolumesTracing.py
Go to the documentation of this file.
1# This file is part of the Peano project. For conditions of distribution and
2# use, please see the copyright notice at www.peano-framework.org
3from peano4.solversteps.ActionSet import ActionSet
4
6
7import jinja2
8
9
12):
13 """!
14
15 Particle tracing over the Finite Volumes solver.
16
17 The action set works for cell-centered Finite Differences, too.
18
19
20 ## Update plug-in points
21
22 We only update cells which are marked with getHasUpdated(). For enclave solvers,
23 we assume that this flag is only set after a cell has updated. Further to that,
24 we also trace in the grid initialisation. If we did not trace there, the first
25 dump prior to the first time step always would yield garbage. It would not even
26 help to add tracing to the plot, as the tracing typically is cell-based, whereas
27 most plotters plug into touchVertexFirstTime(), i.e. they still would dump
28 the particle data prior to the actual tracing.
29
30 The class supports different ways how the data is projected onto tracer
31 attributes (both which data entries and how they are interpolated), as
32 well as various time integratiors. The only important detail here is
33 that the integration happens patch-wisely, i.e. you cannot access any
34 neighbouring patches or intermediate data from the time stepping calculations.
35
36 """
37
39 self,
40 particle_set,
41 solver,
42 project_on_tracer_properties_kernel,
43 projection_kernel_arguments="""
44 marker,
45 {{PATCH_SIZE}},
46 {{NUMBER_OF_UNKNOWNS}}+{{NUMBER_OF_AUXILIARY_VARIABLES}},
47 fineGridCell{{SOLVER_NAME}}Q.value,
48 *p
49""",
50 ):
51 """!
52
53 Set up the tracing
54
55 Setting up the tracing first and foremost means that we have to couple
56 the tracer particle set with the solver, and then we have to specify which
57 projection kernels and which time integrator we want to use. Basically,
58 we have to provide the C++ kernels which shall be used to move the tracer
59 and to write new data into the tracer.
60
61 For the time stepping, Peano's particle toolbox provides some default
62 integrators. If you use stationary tracers, toolbox::particles::StaticPosition
63 (handed in as string) or None do the job. See comment below.
64
65 For the projection of solution data onto the tracer, the file
66 exahype2/fv/Tracer.h provides a few pre-manufactured routines. Again, the
67 projection kernel here requires a plain function call, and you have to
68 ensure that this function call matches the data you actually wanna dump.
69
70
71 ## Arguments
72
73 particle_set: ParticleSet
74
75 time_stepping_kernel: String
76 Function call which does projection onto the particles and the actual
77 position update. Hand in "toolbox::particles::StaticPosition" or None
78 if you don't want the particles' to move.
79
80 project_on_tracer_properties_kernel: String
81 Name of the routine that projects the Finite Volumes representation
82 onto the tracers' properties. Hand in None if you don't want any field
83 properties to be projected on the tracer. This makes sense if the
84 velocity/trajectory of a tracer already encodes all information.
85 Some popular arguments for this routine are enlisted above.
86
87
88 """
89
90 self.tracerDict = {}
91 # Don't overwrite hte dictionary of hte superclass. We'll need it
92 self.tracerDict["PARTICLE"] = particle_set.particle_model.name
93 self.tracerDict["PARTICLES_CONTAINER"] = particle_set.name
94 self.tracerDict["SOLVER_NAME"] = solver._name
95 self.tracerDict["SOLVER_INSTANCE"] = solver.get_name_of_global_instance()
96 self.tracerDict["PATCH_SIZE"] = solver._patch_size
97 self.tracerDict["NUMBER_OF_UNKNOWNS"] = solver._unknowns
98 self.tracerDict["NUMBER_OF_AUXILIARY_VARIABLES"] = solver._auxiliary_variables
99 self.tracerDict["PROJECTION_KERNEL"] = project_on_tracer_properties_kernel
100 self.tracerDict["PROJECTION_KERNEL_ARGUMENTS"] = jinja2.Template(
101 projection_kernel_arguments
102 ).render(**self.tracerDict)
103
104 cell_compute_kernel = """
105if (
106 fineGridCell{{SOLVER_NAME}}CellLabel.getHasUpdated()
107 and
108 (
109 repositories::{{SOLVER_INSTANCE}}.isLastGridSweepOfTimeStep()
110 or
111 repositories::{{SOLVER_INSTANCE}}.getSolverState()=={{SOLVER_NAME}}::SolverState::GridInitialisation
112 )
113) {
114 for (auto* p: localParticles) {
115 if (
116 p->getMoveState()==globaldata::{{PARTICLE}}::MoveState::NotMovedYet and marker.isContained(p->getX())
117 ) {
118 p->getMoveState()==globaldata::{{PARTICLE}}::MoveState::Moved;
119
120
121 {{PROJECTION_KERNEL}}(
122 {{PROJECTION_KERNEL_ARGUMENTS}}
123 );
124 }
125 }
126}
127"""
128
129 touch_vertex_first_time_compute_kernel = """
130 for (auto* p: assignedParticles) {
131 p->setMoveState(globaldata::{{PARTICLE}}::MoveState::NotMovedYet);
132 }
133"""
134
135 super(FiniteVolumesTracing, self).__init__(
136 particle_set,
137 jinja2.Template(cell_compute_kernel).render(**self.tracerDict),
138 jinja2.Template(touch_vertex_first_time_compute_kernel).render(
139 **self.tracerDict
140 ),
141 )
142
143 # def get_body_of_operation(self,operation_name):
144 # result = peano4.toolbox.particles.UpdateParticle_MultiLevelInteraction_StackOfLists_ContiguousParticles.get_body_of_operation(self,operation_name)
145 # if operation_name==ActionSet.OPERATION_BEGIN_TRAVERSAL:
146 # result = self.__Template_BeginIteration.render(**self.tracerDict)
147 # return result
148
150 return __name__.replace(".py", "").replace(".", "_")
151
153 return False
154
155 def get_includes(self):
156 result = jinja2.Template(
157 """
158#include "tarch/multicore/Lock.h"
159#include "peano4/utils/Loop.h"
160#include "vertexdata/{{PARTICLES_CONTAINER}}.h"
161#include "globaldata/{{PARTICLE}}.h"
162#include "repositories/SolverRepository.h"
163#include "exahype2/fv/Tracer.h"
164"""
165 )
166 return (super(FiniteVolumesTracing, self).get_includes()
167 + "\n"
168 + result.render(**self.tracerDict)
169 )
170
171
172 def get_attributes(self):
173 return (
174 super(FiniteVolumesTracing, self).get_attributes()
175 + """
176double _timeStepSize;
177"""
178 )
Particle tracing over the Finite Volumes solver.
user_should_modify_template(self)
Is the user allowed to modify the output.
get_includes(self)
Return include statements that you need.
__init__(self, particle_set, solver, project_on_tracer_properties_kernel, projection_kernel_arguments="""marker,{{PATCH_SIZE}},{{NUMBER_OF_UNKNOWNS}}+{{NUMBER_OF_AUXILIARY_VARIABLES}},fineGridCell{{SOLVER_NAME}}Q.value,*p""")
Set up the tracing.
get_attributes(self)
Return attributes as copied and pasted into the generated class.
Action set (reactions to events)
Definition ActionSet.py:6