ExaDG
Loading...
Searching...
No Matches
flow_rate_calculator.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_FLOW_RATE_CALCULATOR_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_FLOW_RATE_CALCULATOR_H_
24
25// ExaDG
26#include <exadg/matrix_free/integrators.h>
27#include <exadg/utilities/print_functions.h>
28
29namespace ExaDG
30{
31namespace IncNS
32{
33template<int dim>
34struct FlowRateCalculatorData
35{
36 FlowRateCalculatorData() : calculate(false), write_to_file(false), filename("flow_rate")
37 {
38 }
39
40 void
41 print(dealii::ConditionalOStream & pcout)
42 {
43 if(calculate == true)
44 {
45 pcout << " Flow rate calculator:" << std::endl;
46
47 print_parameter(pcout, "Calculate flow rate", calculate);
48 print_parameter(pcout, "Write results to file", write_to_file);
49 if(write_to_file == true)
50 {
51 print_parameter(pcout, "Directory", directory);
52 print_parameter(pcout, "Filename", filename);
53 }
54 }
55 }
56
57 // calculate?
58 bool calculate;
59
60 // write results to file?
61 bool write_to_file;
62
63 // directory and filename
64 std::string directory;
65 std::string filename;
66};
67
68
69/*
70 * This class calculates a vector of flow rates where the different entries of the vector correspond
71 * to different boundary IDs, i.e., one outflow boundary may only consist of faces with the same
72 * boundary ID. This class is intended to be used in cases where the number of outflow boundaries is
73 * very large with a typical use case being the human lung (where the number of outflow boundaries
74 * is 2^{N-1} with N being the number of airway generations). The reason behind is that calculating
75 * the flow rate requires global communication since different processors may share the same outflow
76 * boundary and the implementation in this class only requires one communication for all outflow
77 * boundaries.
78 *
79 * Note: For a small number of outflow boundaries (for example 1-10), the more modular class
80 * MeanVelocityCalculator should be used rather than this specialized implementation.
81 */
82template<int dim, typename Number>
83class FlowRateCalculator
84{
85public:
86 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
87
88 typedef CellIntegrator<dim, dim, Number> CellIntegratorU;
89 typedef FaceIntegrator<dim, dim, Number> FaceIntegratorU;
90
91 typedef dealii::VectorizedArray<Number> scalar;
92
93 FlowRateCalculator(dealii::MatrixFree<dim, Number> const & matrix_free_in,
94 unsigned int const dof_index_in,
95 unsigned int const quad_index_in,
96 FlowRateCalculatorData<dim> const & data_in,
97 MPI_Comm const & mpi_comm_in);
98
99 Number
100 calculate_flow_rates(VectorType const & velocity,
101 double const & time,
102 std::map<dealii::types::boundary_id, Number> & flow_rates);
103
104
105private:
106 void
107 write_output(Number const & value, double const & time, std::string const & name);
108
109 void
110 do_calculate_flow_rates(VectorType const & velocity,
111 std::map<dealii::types::boundary_id, Number> & flow_rates);
112
113 FlowRateCalculatorData<dim> const & data;
114 dealii::MatrixFree<dim, Number> const & matrix_free;
115 unsigned int dof_index, quad_index;
116 bool clear_files;
117
118 MPI_Comm const mpi_comm;
119};
120
121} // namespace IncNS
122} // namespace ExaDG
123
124#endif /* EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_FLOW_RATE_CALCULATOR_H_ */
Definition driver.cpp:33
Definition flow_rate_calculator.h:35