186 * We declare
class that describes the boundary conditions and
initial one for velocity:
193 * class Velocity:
public Function<dim> {
195 * Velocity(
const double initial_time = 0.0);
198 *
const unsigned int component = 0)
const override;
206 * Velocity<dim>::Velocity(
const double initial_time):
Function<dim>(dim, initial_time) {}
210 *
double Velocity<dim>::value(
const Point<dim>& p,
const unsigned int component)
const {
212 *
if(component == 0) {
213 *
const double Um = 1.5;
214 *
const double H = 4.1;
216 *
return 4.0*Um*p(1)*(H - p(1))/(H*H);
227 *
for(
unsigned int i = 0; i < dim; ++i)
228 * values[i] =
value(p, i);
234 * We
do the same
for the pressure
241 *
class Pressure:
public Function<dim> {
243 * Pressure(
const double initial_time = 0.0);
246 *
const unsigned int component = 0)
const override;
251 * Pressure<dim>::Pressure(
const double initial_time):
Function<dim>(1, initial_time) {}
255 *
double Pressure<dim>::value(
const Point<dim>& p,
const unsigned int component)
const {
259 *
return 22.0 - p(0);
266<a name=
"ann-navier_stokes_TRBDF2_DG.cc"></a>
267<h1>Annotated version of navier_stokes_TRBDF2_DG.cc</h1>
277 * We start by including all the necessary deal.II header files and some
C++
284 * #include <deal.II/base/quadrature_lib.h>
285 * #include <deal.II/base/multithread_info.h>
286 * #include <deal.II/base/thread_management.h>
287 * #include <deal.II/base/work_stream.h>
288 * #include <deal.II/base/
parallel.h>
289 * #include <deal.II/base/utilities.h>
290 * #include <deal.II/base/conditional_ostream.h>
292 * #include <deal.II/lac/vector.h>
293 * #include <deal.II/lac/solver_cg.h>
294 * #include <deal.II/lac/precondition.h>
295 * #include <deal.II/lac/solver_gmres.h>
296 * #include <deal.II/lac/affine_constraints.h>
298 * #include <deal.II/grid/
tria.h>
299 * #include <deal.II/grid/grid_generator.h>
300 * #include <deal.II/grid/grid_tools.h>
301 * #include <deal.II/grid/grid_refinement.h>
302 * #include <deal.II/grid/tria_accessor.h>
303 * #include <deal.II/grid/tria_iterator.h>
304 * #include <deal.II/distributed/grid_refinement.h>
306 * #include <deal.II/dofs/dof_handler.h>
307 * #include <deal.II/dofs/dof_accessor.h>
308 * #include <deal.II/dofs/dof_tools.h>
310 * #include <deal.II/fe/fe_q.h>
311 * #include <deal.II/fe/fe_dgq.h>
312 * #include <deal.II/fe/fe_values.h>
313 * #include <deal.II/fe/fe_tools.h>
314 * #include <deal.II/fe/fe_system.h>
316 * #include <deal.II/numerics/matrix_tools.h>
317 * #include <deal.II/numerics/vector_tools.h>
318 * #include <deal.II/numerics/data_out.h>
322 * #include <iostream>
324 * #include <deal.II/matrix_free/matrix_free.h>
325 * #include <deal.II/matrix_free/operators.h>
326 * #include <deal.II/matrix_free/fe_evaluation.h>
327 * #include <deal.II/fe/component_mask.h>
329 * #include <deal.II/base/timer.h>
330 * #include <deal.II/distributed/solution_transfer.h>
331 * #include <deal.II/numerics/error_estimator.h>
333 * #include <deal.II/multigrid/multigrid.h>
334 * #include <deal.II/multigrid/mg_transfer_matrix_free.h>
335 * #include <deal.II/multigrid/mg_tools.h>
336 * #include <deal.II/multigrid/mg_coarse.h>
337 * #include <deal.II/multigrid/mg_smoother.h>
338 * #include <deal.II/multigrid/mg_matrix.h>
340 * #include <deal.II/meshworker/
mesh_loop.h>
342 * #include
"runtime_parameters.h"
343 * #include
"equation_data.h"
348 *
template<
int dim,
typename Number,
typename VectorizedArrayType>
353 *
const unsigned int&,
354 *
const std::pair<unsigned int, unsigned int>&)>& cell_operation,
357 *
const unsigned int&,
358 *
const std::pair<unsigned int, unsigned int>&)>& face_operation,
361 *
const unsigned int&,
362 *
const std::pair<unsigned int, unsigned int>&)>& boundary_operation,
363 *
const unsigned int dof_no = 0) {
371 *
const unsigned int dummy = 0;
373 * matrix_free.
loop(cell_operation, face_operation, boundary_operation,
374 * diagonal_global, dummy,
false,
382 * We include the code in a suitable
namespace:
388 *
namespace NS_TRBDF2 {
393 * The following
class is an auxiliary one for post-processing of the vorticity
403 * std::vector<
Vector<double>>& computed_quantities)
const override;
405 *
virtual std::vector<std::string> get_names()
const override;
407 *
virtual std::vector<DataComponentInterpretation::DataComponentInterpretation>
408 * get_data_component_interpretation()
const override;
410 *
virtual UpdateFlags get_needed_update_flags()
const override;
415 * This function evaluates the vorticty in both 2D and 3D cases
424 *
const unsigned int n_quadrature_points = inputs.
solution_values.size();
428 *
Assert(computed_quantities.size() == n_quadrature_points, ExcInternalError());
433 *
Assert(computed_quantities[0].size() == 1, ExcInternalError());
436 *
Assert(computed_quantities[0].size() == dim, ExcInternalError());
441 *
for(
unsigned int q = 0; q < n_quadrature_points; ++q)
445 *
for(
unsigned int q = 0; q < n_quadrature_points; ++q) {
455 * This auxiliary function is required by the base
class DataProcessor and simply
456 * sets the name
for the output file
463 * std::vector<std::string> PostprocessorVorticity<dim>::get_names()
const {
464 * std::vector<std::string> names;
465 * names.emplace_back(
"vorticity");
467 * names.emplace_back(
"vorticity");
468 * names.emplace_back(
"vorticity");
476 * This auxiliary function is required by the base
class DataProcessor and simply
477 * specifies
if the vorticity is a
scalar (2D) or a vector (3D)
484 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
485 * PostprocessorVorticity<dim>::get_data_component_interpretation()
const {
486 * std::vector<DataComponentInterpretation::DataComponentInterpretation> interpretation;
495 *
return interpretation;
500 * This auxiliary function is required by the base
class DataProcessor and simply
501 * sets which variables have to updated (only the gradients)
508 *
UpdateFlags PostprocessorVorticity<dim>::get_needed_update_flags()
const {
515 * The following structs are auxiliary objects
for mesh refinement. ScratchData simply sets
523 *
struct ScratchData {
525 *
const unsigned int quadrature_degree,
526 *
const UpdateFlags update_flags): fe_values(fe,
QGauss<dim>(quadrature_degree), update_flags) {}
528 * ScratchData(
const ScratchData<dim>& scratch_data): fe_values(scratch_data.fe_values.get_fe(),
529 * scratch_data.fe_values.get_quadrature(),
530 * scratch_data.fe_values.get_update_flags()) {}
537 * CopyData simply sets the cell
index
546 * CopyData(
const CopyData &) =
default;
557 * @sect{ <code>NavierStokesProjectionOperator::NavierStokesProjectionOperator</code> }
561 * The following
class sets effecively the weak formulation of the problems for the different stages
562 * and
for both velocity and pressure.
563 * The
template parameters are the dimnesion of the problem, the polynomial degree
for the pressure,
564 * the polynomial degree
for the velocity, the number of quadrature points
for integrals
for the pressure step,
565 * the number of quadrature points
for integrals
for the velocity step, the type of vector
for storage and the type
566 * of floating
point data (in general
double or
float for preconditioners structures
if desired).
572 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
575 *
using Number =
typename Vec::value_type;
577 * NavierStokesProjectionOperator();
579 * NavierStokesProjectionOperator(RunTimeParameters::Data_Storage& data);
581 *
void set_dt(
const double time_step);
583 *
void set_TR_BDF2_stage(
const unsigned int stage);
585 *
void set_NS_stage(
const unsigned int stage);
587 *
void set_u_extr(
const Vec& src);
589 *
void vmult_rhs_velocity(Vec& dst,
const std::vector<Vec>& src)
const;
591 *
void vmult_rhs_pressure(Vec& dst,
const std::vector<Vec>& src)
const;
593 *
void vmult_grad_p_projection(Vec& dst,
const Vec& src)
const;
607 *
unsigned int TR_BDF2_stage;
608 *
unsigned int NS_stage;
611 *
virtual void apply_add(Vec& dst,
const Vec& src)
const override;
616 *
const double a21 = 0.5;
617 *
const double a22 = 0.5;
620 *
const double theta_v = 1.0;
621 *
const double theta_p = 1.0;
622 *
const double C_p = 1.0*(fe_degree_p + 1)*(fe_degree_p + 1);
623 *
const double C_u = 1.0*(fe_degree_v + 1)*(fe_degree_v + 1);
627 * EquationData::Velocity<dim> vel_boundary_inflow;
633 *
const std::vector<Vec>& src,
634 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
637 *
const std::vector<Vec>& src,
638 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
641 *
const std::vector<Vec>& src,
642 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
646 *
const std::vector<Vec>& src,
647 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
650 *
const std::vector<Vec>& src,
651 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
654 *
const std::vector<Vec>& src,
655 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
660 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
664 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
668 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
673 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
677 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
681 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
686 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
690 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
694 *
const unsigned int& src,
695 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
698 *
const unsigned int& src,
699 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
702 *
const unsigned int& src,
703 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
707 *
const unsigned int& src,
708 *
const std::pair<unsigned int, unsigned int>& cell_range)
const;
711 *
const unsigned int& src,
712 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
715 *
const unsigned int& src,
716 *
const std::pair<unsigned int, unsigned int>& face_range)
const;
722 * We start with the
default constructor. It is important
for MultiGrid, so it is fundamental
723 * to properly
set the parameters of the time scheme.
729 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
730 * NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
731 * NavierStokesProjectionOperator():
733 * a32(a31), a33(1.0/(2.0 -
gamma)), TR_BDF2_stage(1), NS_stage(1), u_extr() {}
738 * We focus now on the constructor with runtime parameters storage
744 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
745 * NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
746 * NavierStokesProjectionOperator(RunTimeParameters::Data_Storage& data):
749 * a32(a31), a33(1.0/(2.0 -
gamma)), TR_BDF2_stage(1), NS_stage(1), u_extr(),
750 * vel_boundary_inflow(data.initial_time) {}
755 * Setter of time-step (called by
Multigrid and in
case a smaller time-step towards the end is needed)
761 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
762 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
763 * set_dt(
const double time_step) {
770 * Setter of TR-BDF2 stage (
this can be known only during the effective execution
771 * and so it has to be demanded to the
class that really solves the problem)
777 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
778 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
779 * set_TR_BDF2_stage(
const unsigned int stage) {
781 *
Assert(stage > 0, ExcInternalError());
783 * TR_BDF2_stage = stage;
789 * Setter of NS stage (
this can be known only during the effective execution
790 * and so it has to be demanded to the
class that really solves the problem)
796 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
797 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
798 * set_NS_stage(
const unsigned int stage) {
800 *
Assert(stage > 0, ExcInternalError());
808 * Setter of extrapolated velocity
for different stages
814 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
815 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
816 * set_u_extr(
const Vec& src) {
818 * u_extr.update_ghost_values();
824 * We are in a DG-
MatrixFree framework, so it is convenient to compute separately cell contribution,
825 *
internal faces contributions and boundary faces contributions. We start by
826 * assembling the rhs cell term
for the velocity.
832 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
833 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
836 *
const std::vector<Vec>& src,
837 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
838 *
if(TR_BDF2_stage == 1) {
845 * phi_old_extr(data, 0);
849 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
853 * phi_old.reinit(cell);
858 * phi_old_extr.reinit(cell);
860 * phi_old_press.reinit(cell);
865 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
866 *
const auto& u_n = phi_old.get_value(q);
867 *
const auto& grad_u_n = phi_old.get_gradient(q);
868 *
const auto& u_n_gamma_ov_2 = phi_old_extr.get_value(q);
869 *
const auto& tensor_product_u_n =
outer_product(u_n, u_n_gamma_ov_2);
870 *
const auto& p_n = phi_old_press.get_value(q);
871 *
auto p_n_times_identity = tensor_product_u_n;
872 * p_n_times_identity = 0;
873 *
for(
unsigned int d = 0;
d < dim; ++
d)
874 * p_n_times_identity[d][d] = p_n;
876 * phi.submit_value(1.0/(gamma*dt)*u_n, q);
878 * phi.submit_gradient(-a21/Re*grad_u_n + a21*tensor_product_u_n + p_n_times_identity, q);
894 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
895 * phi_old.reinit(cell);
897 * phi_int.reinit(cell);
899 * phi_old_press.reinit(cell);
904 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
905 *
const auto& u_n = phi_old.get_value(q);
906 *
const auto& grad_u_n = phi_old.get_gradient(q);
907 *
const auto& u_n_gamma = phi_int.get_value(q);
908 *
const auto& grad_u_n_gamma = phi_int.get_gradient(q);
910 *
const auto& tensor_product_u_n_gamma =
outer_product(u_n_gamma, u_n_gamma);
911 *
const auto& p_n = phi_old_press.get_value(q);
912 *
auto p_n_times_identity = tensor_product_u_n;
913 * p_n_times_identity = 0;
914 *
for(
unsigned int d = 0;
d < dim; ++
d)
915 * p_n_times_identity[d][d] = p_n;
917 * phi.submit_value(1.0/((1.0 - gamma)*dt)*u_n_gamma, q);
918 * phi.submit_gradient(a32*tensor_product_u_n_gamma + a31*tensor_product_u_n -
919 * a32/Re*grad_u_n_gamma - a31/Re*grad_u_n + p_n_times_identity, q);
929 * The followinf function assembles rhs face term
for the velocity
935 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
936 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
939 *
const std::vector<Vec>& src,
940 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
941 *
if(TR_BDF2_stage == 1) {
946 * phi_m(data,
false, 0),
947 * phi_old_p(data,
true, 0),
948 * phi_old_m(data,
false, 0),
949 * phi_old_extr_p(data,
true, 0),
950 * phi_old_extr_m(data,
false, 0);
952 * phi_old_press_m(data,
false, 1);
955 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
956 * phi_old_p.reinit(face);
958 * phi_old_m.reinit(face);
960 * phi_old_extr_p.reinit(face);
962 * phi_old_extr_m.reinit(face);
964 * phi_old_press_p.reinit(face);
966 * phi_old_press_m.reinit(face);
968 * phi_p.reinit(face);
969 * phi_m.reinit(face);
972 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
973 *
const auto& n_plus = phi_p.get_normal_vector(q);
977 *
const auto& avg_grad_u_old = 0.5*(phi_old_p.get_gradient(q) + phi_old_m.get_gradient(q));
978 *
const auto& avg_tensor_product_u_n = 0.5*(
outer_product(phi_old_p.get_value(q), phi_old_extr_p.get_value(q)) +
979 *
outer_product(phi_old_m.get_value(q), phi_old_extr_m.get_value(q)));
980 *
const auto& avg_p_old = 0.5*(phi_old_press_p.get_value(q) + phi_old_press_m.get_value(q));
982 * phi_p.submit_value((a21/Re*avg_grad_u_old - a21*avg_tensor_product_u_n)*n_plus - avg_p_old*n_plus, q);
983 * phi_m.submit_value(-(a21/Re*avg_grad_u_old - a21*avg_tensor_product_u_n)*n_plus + avg_p_old*n_plus, q);
992 * phi_m(data,
false, 0),
993 * phi_old_p(data,
true, 0),
994 * phi_old_m(data,
false, 0),
995 * phi_int_p(data,
true, 0),
996 * phi_int_m(data,
false, 0);
998 * phi_old_press_m(data,
false, 1);
1001 *
for(
unsigned int face = face_range.first; face < face_range.second; ++ face) {
1002 * phi_old_p.reinit(face);
1004 * phi_old_m.reinit(face);
1006 * phi_int_p.reinit(face);
1008 * phi_int_m.reinit(face);
1010 * phi_old_press_p.reinit(face);
1012 * phi_old_press_m.reinit(face);
1014 * phi_p.reinit(face);
1015 * phi_m.reinit(face);
1018 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
1019 *
const auto& n_plus = phi_p.get_normal_vector(q);
1021 *
const auto& avg_grad_u_old = 0.5*(phi_old_p.get_gradient(q) + phi_old_m.get_gradient(q));
1022 *
const auto& avg_grad_u_int = 0.5*(phi_int_p.get_gradient(q) + phi_int_m.get_gradient(q));
1023 *
const auto& avg_tensor_product_u_n = 0.5*(
outer_product(phi_old_p.get_value(q), phi_old_p.get_value(q)) +
1024 *
outer_product(phi_old_m.get_value(q), phi_old_m.get_value(q)));
1025 *
const auto& avg_tensor_product_u_n_gamma = 0.5*(
outer_product(phi_int_p.get_value(q), phi_int_p.get_value(q)) +
1026 *
outer_product(phi_int_m.get_value(q), phi_int_m.get_value(q)));
1027 *
const auto& avg_p_old = 0.5*(phi_old_press_p.get_value(q) + phi_old_press_m.get_value(q));
1029 * phi_p.submit_value((a31/Re*avg_grad_u_old + a32/Re*avg_grad_u_int -
1030 * a31*avg_tensor_product_u_n - a32*avg_tensor_product_u_n_gamma)*n_plus - avg_p_old*n_plus, q);
1031 * phi_m.submit_value(-(a31/Re*avg_grad_u_old + a32/Re*avg_grad_u_int -
1032 * a31*avg_tensor_product_u_n - a32*avg_tensor_product_u_n_gamma)*n_plus + avg_p_old*n_plus, q);
1043 * The followinf function assembles rhs boundary term
for the velocity
1049 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1050 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1053 *
const std::vector<Vec>& src,
1054 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1055 *
if(TR_BDF2_stage == 1) {
1059 * phi_old(data,
true, 0),
1060 * phi_old_extr(data,
true, 0);
1064 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1065 * phi_old.reinit(face);
1067 * phi_old_extr.reinit(face);
1069 * phi_old_press.reinit(face);
1074 *
const auto coef_jump = (
boundary_id == 1) ? 0.0 : C_u*
std::
abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
1075 *
const double aux_coeff = (
boundary_id == 1) ? 0.0 : 1.0;
1078 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1079 *
const auto& n_plus = phi.get_normal_vector(q);
1081 *
const auto& grad_u_old = phi_old.get_gradient(q);
1082 *
const auto& tensor_product_u_n =
outer_product(phi_old.get_value(q), phi_old_extr.get_value(q));
1083 *
const auto& p_old = phi_old_press.get_value(q);
1084 *
const auto& point_vectorized = phi.quadrature_point(q);
1086 *
if(boundary_id == 0) {
1087 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
1091 *
for(
unsigned int d = 0;
d < dim; ++
d)
1092 * point[d] = point_vectorized[d][v];
1093 *
for(
unsigned int d = 0;
d < dim; ++
d)
1094 * u_int_m[d][v] = vel_boundary_inflow.value(point, d);
1097 *
const auto tensor_product_u_int_m =
outer_product(u_int_m, phi_old_extr.get_value(q));
1100 * phi.submit_value((a21/Re*grad_u_old - a21*tensor_product_u_n)*n_plus - p_old*n_plus +
1101 * a22/Re*2.0*coef_jump*u_int_m -
1102 * aux_coeff*a22*tensor_product_u_int_m*n_plus + a22*lambda*u_int_m, q);
1103 * phi.submit_normal_derivative(-aux_coeff*theta_v*a22/Re*u_int_m, q);
1112 * phi_old(data,
true, 0),
1113 * phi_int(data,
true, 0),
1114 * phi_int_extr(data,
true, 0);
1118 *
for(
unsigned int face = face_range.first; face < face_range.second; ++ face) {
1119 * phi_old.reinit(face);
1121 * phi_int.reinit(face);
1123 * phi_old_press.reinit(face);
1125 * phi_int_extr.reinit(face);
1130 *
const auto coef_jump = (
boundary_id == 1) ? 0.0 : C_u*
std::
abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
1131 *
const double aux_coeff = (
boundary_id == 1) ? 0.0 : 1.0;
1134 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1135 *
const auto& n_plus = phi.get_normal_vector(q);
1137 *
const auto& grad_u_old = phi_old.get_gradient(q);
1138 *
const auto& grad_u_int = phi_int.get_gradient(q);
1139 *
const auto& tensor_product_u_n =
outer_product(phi_old.get_value(q), phi_old.get_value(q));
1140 *
const auto& tensor_product_u_n_gamma =
outer_product(phi_int.get_value(q), phi_int.get_value(q));
1141 *
const auto& p_old = phi_old_press.get_value(q);
1142 *
const auto& point_vectorized = phi.quadrature_point(q);
1144 *
if(boundary_id == 0) {
1145 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
1147 *
for(
unsigned int d = 0;
d < dim; ++
d)
1148 * point[d] = point_vectorized[d][v];
1149 *
for(
unsigned int d = 0;
d < dim; ++
d)
1150 * u_m[d][v] = vel_boundary_inflow.value(point, d);
1153 *
const auto tensor_product_u_m =
outer_product(u_m, phi_int_extr.get_value(q));
1156 * phi.submit_value((a31/Re*grad_u_old + a32/Re*grad_u_int -
1157 * a31*tensor_product_u_n - a32*tensor_product_u_n_gamma)*n_plus - p_old*n_plus +
1158 * a33/Re*2.0*coef_jump*u_m -
1159 * aux_coeff*a33*tensor_product_u_m*n_plus + a33*lambda*u_m, q);
1160 * phi.submit_normal_derivative(-aux_coeff*theta_v*a33/Re*u_m, q);
1170 * Put together all the previous steps
for velocity. This is done automatically by the
loop function of
'MatrixFree' class
1176 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1177 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1178 * vmult_rhs_velocity(Vec& dst,
const std::vector<Vec>& src)
const {
1179 *
for(
auto& vec : src)
1180 * vec.update_ghost_values();
1182 * this->data->
loop(&NavierStokesProjectionOperator::assemble_rhs_cell_term_velocity,
1183 * &NavierStokesProjectionOperator::assemble_rhs_face_term_velocity,
1184 * &NavierStokesProjectionOperator::assemble_rhs_boundary_term_velocity,
1185 *
this, dst, src,
true,
1193 * Now we focus on computing the rhs
for the projection step
for the pressure with the same ratio.
1194 * The following function assembles rhs cell term
for the pressure
1200 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1201 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1204 *
const std::vector<Vec>& src,
1205 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1209 * phi_old(data, 1, 1);
1212 *
const double coeff = (TR_BDF2_stage == 1) ? 1.0e6*gamma*dt*gamma*dt : 1.0e6*(1.0 -
gamma)*dt*(1.0 -
gamma)*dt;
1214 *
const double coeff_2 = (TR_BDF2_stage == 1) ? gamma*dt : (1.0 -
gamma)*dt;
1217 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1218 * phi_proj.reinit(cell);
1220 * phi_old.reinit(cell);
1225 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1226 *
const auto& u_star_star = phi_proj.get_value(q);
1227 *
const auto& p_old = phi_old.get_value(q);
1229 * phi.submit_value(1.0/coeff*p_old, q);
1230 * phi.submit_gradient(1.0/coeff_2*u_star_star, q);
1239 * The following function assembles rhs face term
for the pressure
1245 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1246 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1249 *
const std::vector<Vec>& src,
1250 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1253 * phi_m(data,
false, 1, 1);
1255 * phi_proj_m(data,
false, 0, 1);
1257 *
const double coeff = (TR_BDF2_stage == 1) ? 1.0/(gamma*dt) : 1.0/((1.0 -
gamma)*dt);
1260 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1261 * phi_proj_p.reinit(face);
1263 * phi_proj_m.reinit(face);
1265 * phi_p.reinit(face);
1266 * phi_m.reinit(face);
1269 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
1270 *
const auto& n_plus = phi_p.get_normal_vector(q);
1271 *
const auto& avg_u_star_star = 0.5*(phi_proj_p.get_value(q) + phi_proj_m.get_value(q));
1273 * phi_p.submit_value(-coeff*scalar_product(avg_u_star_star, n_plus), q);
1274 * phi_m.submit_value(coeff*scalar_product(avg_u_star_star, n_plus), q);
1284 * The following function assembles rhs boundary term
for the pressure
1290 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1291 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1294 *
const std::vector<Vec>& src,
1295 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1300 *
const double coeff = (TR_BDF2_stage == 1) ? 1.0/(gamma*dt) : 1.0/((1.0 -
gamma)*dt);
1303 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1304 * phi_proj.reinit(face);
1309 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1310 *
const auto& n_plus = phi.get_normal_vector(q);
1312 * phi.submit_value(-coeff*scalar_product(phi_proj.get_value(q), n_plus), q);
1321 * Put together all the previous steps
for pressure
1327 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1328 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1329 * vmult_rhs_pressure(Vec& dst,
const std::vector<Vec>& src)
const {
1330 *
for(
auto& vec : src)
1331 * vec.update_ghost_values();
1333 * this->data->
loop(&NavierStokesProjectionOperator::assemble_rhs_cell_term_pressure,
1334 * &NavierStokesProjectionOperator::assemble_rhs_face_term_pressure,
1335 * &NavierStokesProjectionOperator::assemble_rhs_boundary_term_pressure,
1336 *
this, dst, src,
true,
1344 * Now we need to build the
'matrices', i.e. the bilinear forms. We start by
1345 * assembling the cell term
for the velocity
1351 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1352 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1356 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1357 *
if(TR_BDF2_stage == 1) {
1361 * phi_old_extr(data, 0);
1364 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1367 * phi_old_extr.reinit(cell);
1371 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1372 *
const auto& u_int = phi.get_value(q);
1373 *
const auto& grad_u_int = phi.get_gradient(q);
1374 *
const auto& u_n_gamma_ov_2 = phi_old_extr.get_value(q);
1375 *
const auto& tensor_product_u_int =
outer_product(u_int, u_n_gamma_ov_2);
1377 * phi.submit_value(1.0/(gamma*dt)*u_int, q);
1378 * phi.submit_gradient(-a22*tensor_product_u_int + a22/Re*grad_u_int, q);
1386 * phi_int_extr(data, 0);
1389 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1392 * phi_int_extr.reinit(cell);
1396 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1397 *
const auto& u_curr = phi.get_value(q);
1398 *
const auto& grad_u_curr = phi.get_gradient(q);
1399 *
const auto& u_n1_int = phi_int_extr.get_value(q);
1400 *
const auto& tensor_product_u_curr =
outer_product(u_curr, u_n1_int);
1402 * phi.submit_value(1.0/((1.0 - gamma)*dt)*u_curr, q);
1403 * phi.submit_gradient(-a33*tensor_product_u_curr + a33/Re*grad_u_curr, q);
1413 * The following function assembles face term
for the velocity
1419 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1420 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1424 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1425 *
if(TR_BDF2_stage == 1) {
1428 * phi_m(data,
false, 0),
1429 * phi_old_extr_p(data,
true, 0),
1430 * phi_old_extr_m(data,
false, 0);
1433 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1434 * phi_p.reinit(face);
1436 * phi_m.reinit(face);
1438 * phi_old_extr_p.reinit(face);
1440 * phi_old_extr_m.reinit(face);
1443 *
const auto coef_jump = C_u*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
1444 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
1447 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
1448 *
const auto& n_plus = phi_p.get_normal_vector(q);
1450 *
const auto& avg_grad_u_int = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
1451 *
const auto& jump_u_int = phi_p.get_value(q) - phi_m.get_value(q);
1452 *
const auto& avg_tensor_product_u_int = 0.5*(
outer_product(phi_p.get_value(q), phi_old_extr_p.get_value(q)) +
1453 *
outer_product(phi_m.get_value(q), phi_old_extr_m.get_value(q)));
1455 *
std::abs(scalar_product(phi_old_extr_m.get_value(q), n_plus)));
1457 * phi_p.submit_value(a22/Re*(-avg_grad_u_int*n_plus + coef_jump*jump_u_int) +
1458 * a22*avg_tensor_product_u_int*n_plus + 0.5*a22*lambda*jump_u_int, q);
1459 * phi_m.submit_value(-a22/Re*(-avg_grad_u_int*n_plus + coef_jump*jump_u_int) -
1460 * a22*avg_tensor_product_u_int*n_plus - 0.5*a22*lambda*jump_u_int, q);
1461 * phi_p.submit_normal_derivative(-theta_v*a22/Re*0.5*jump_u_int, q);
1462 * phi_m.submit_normal_derivative(-theta_v*a22/Re*0.5*jump_u_int, q);
1471 * phi_m(data,
false, 0),
1472 * phi_extr_p(data,
true, 0),
1473 * phi_extr_m(data,
false, 0);
1476 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1477 * phi_p.reinit(face);
1479 * phi_m.reinit(face);
1481 * phi_extr_p.reinit(face);
1483 * phi_extr_m.reinit(face);
1486 *
const auto coef_jump = C_u*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
1487 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
1490 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
1491 *
const auto& n_plus = phi_p.get_normal_vector(q);
1493 *
const auto& avg_grad_u = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
1494 *
const auto& jump_u = phi_p.get_value(q) - phi_m.get_value(q);
1495 *
const auto& avg_tensor_product_u = 0.5*(
outer_product(phi_p.get_value(q), phi_extr_p.get_value(q)) +
1496 *
outer_product(phi_m.get_value(q), phi_extr_m.get_value(q)));
1498 *
std::abs(scalar_product(phi_extr_m.get_value(q), n_plus)));
1500 * phi_p.submit_value(a33/Re*(-avg_grad_u*n_plus + coef_jump*jump_u) +
1501 * a33*avg_tensor_product_u*n_plus + 0.5*a33*lambda*jump_u, q);
1502 * phi_m.submit_value(-a33/Re*(-avg_grad_u*n_plus + coef_jump*jump_u) -
1503 * a33*avg_tensor_product_u*n_plus - 0.5*a33*lambda*jump_u, q);
1504 * phi_p.submit_normal_derivative(-theta_v*a33/Re*0.5*jump_u, q);
1505 * phi_m.submit_normal_derivative(-theta_v*a33/Re*0.5*jump_u, q);
1516 * The following function assembles boundary term
for the velocity
1522 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1523 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1527 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1528 *
if(TR_BDF2_stage == 1) {
1531 * phi_old_extr(data,
true, 0);
1534 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1537 * phi_old_extr.reinit(face);
1541 *
const auto coef_jump = C_u*
std::abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
1545 *
if(boundary_id != 1) {
1546 *
const double coef_trasp = 0.0;
1549 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1550 *
const auto& n_plus = phi.get_normal_vector(q);
1551 *
const auto& grad_u_int = phi.get_gradient(q);
1552 *
const auto& u_int = phi.get_value(q);
1553 *
const auto& tensor_product_u_int =
outer_product(phi.get_value(q), phi_old_extr.get_value(q));
1554 *
const auto&
lambda =
std::abs(scalar_product(phi_old_extr.get_value(q), n_plus));
1556 * phi.submit_value(a22/Re*(-grad_u_int*n_plus + 2.0*coef_jump*u_int) +
1557 * a22*coef_trasp*tensor_product_u_int*n_plus + a22*lambda*u_int, q);
1558 * phi.submit_normal_derivative(-theta_v*a22/Re*u_int, q);
1564 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1565 *
const auto& n_plus = phi.get_normal_vector(q);
1566 *
const auto& grad_u_int = phi.get_gradient(q);
1567 *
const auto& u_int = phi.get_value(q);
1568 *
const auto&
lambda =
std::abs(scalar_product(phi_old_extr.get_value(q), n_plus));
1570 *
const auto& point_vectorized = phi.quadrature_point(q);
1571 *
auto u_int_m = u_int;
1572 *
auto grad_u_int_m = grad_u_int;
1573 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
1575 *
for(
unsigned int d = 0;
d < dim; ++
d)
1576 * point[d] = point_vectorized[d][v];
1578 * u_int_m[1][v] = -u_int_m[1][v];
1580 * grad_u_int_m[0][0][v] = -grad_u_int_m[0][0][v];
1581 * grad_u_int_m[0][1][v] = -grad_u_int_m[0][1][v];
1584 * phi.submit_value(a22/Re*(-(0.5*(grad_u_int + grad_u_int_m))*n_plus + coef_jump*(u_int - u_int_m)) +
1585 * a22*
outer_product(0.5*(u_int + u_int_m), phi_old_extr.get_value(q))*n_plus +
1586 * a22*0.5*
lambda*(u_int - u_int_m), q);
1587 * phi.submit_normal_derivative(-theta_v*a22/Re*(u_int - u_int_m), q);
1596 * phi_extr(data,
true, 0);
1599 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1602 * phi_extr.reinit(face);
1606 *
const auto coef_jump = C_u*
std::abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
1608 *
if(boundary_id != 1) {
1609 *
const double coef_trasp = 0.0;
1612 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1613 *
const auto& n_plus = phi.get_normal_vector(q);
1614 *
const auto& grad_u = phi.get_gradient(q);
1615 *
const auto& u = phi.get_value(q);
1616 *
const auto& tensor_product_u =
outer_product(phi.get_value(q), phi_extr.get_value(q));
1617 *
const auto&
lambda =
std::abs(scalar_product(phi_extr.get_value(q), n_plus));
1619 * phi.submit_value(a33/Re*(-grad_u*n_plus + 2.0*coef_jump*u) +
1620 * a33*coef_trasp*tensor_product_u*n_plus + a33*lambda*u, q);
1621 * phi.submit_normal_derivative(-theta_v*a33/Re*u, q);
1627 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1628 *
const auto& n_plus = phi.get_normal_vector(q);
1629 *
const auto& grad_u = phi.get_gradient(q);
1630 *
const auto& u = phi.get_value(q);
1631 *
const auto&
lambda =
std::abs(scalar_product(phi_extr.get_value(q), n_plus));
1633 *
const auto& point_vectorized = phi.quadrature_point(q);
1635 *
auto grad_u_m = grad_u;
1636 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
1638 *
for(
unsigned int d = 0;
d < dim; ++
d)
1639 * point[d] = point_vectorized[d][v];
1641 * u_m[1][v] = -u_m[1][v];
1643 * grad_u_m[0][0][v] = -grad_u_m[0][0][v];
1644 * grad_u_m[0][1][v] = -grad_u_m[0][1][v];
1647 * phi.submit_value(a33/Re*(-(0.5*(grad_u + grad_u_m))*n_plus + coef_jump*(u - u_m)) +
1648 * a33*
outer_product(0.5*(u + u_m), phi_extr.get_value(q))*n_plus + a33*0.5*
lambda*(u - u_m), q);
1649 * phi.submit_normal_derivative(-theta_v*a33/Re*(u - u_m), q);
1660 * Next, we focus on
'matrices' to compute the pressure. We
first assemble cell term
for the pressure
1666 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1667 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1671 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1675 *
const double coeff = (TR_BDF2_stage == 1) ? 1.0e6*gamma*dt*gamma*dt : 1.0e6*(1.0 -
gamma)*dt*(1.0 -
gamma)*dt;
1678 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1683 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1684 * phi.submit_gradient(phi.get_gradient(q), q);
1685 * phi.submit_value(1.0/coeff*phi.get_value(q), q);
1695 * The following function assembles face term
for the pressure
1701 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1702 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1706 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1709 * phi_m(data,
false, 1, 1);
1712 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1713 * phi_p.reinit(face);
1715 * phi_m.reinit(face);
1718 *
const auto coef_jump = C_p*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
1719 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
1722 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
1723 *
const auto& n_plus = phi_p.get_normal_vector(q);
1725 *
const auto& avg_grad_pres = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
1726 *
const auto& jump_pres = phi_p.get_value(q) - phi_m.get_value(q);
1728 * phi_p.submit_value(-scalar_product(avg_grad_pres, n_plus) + coef_jump*jump_pres, q);
1729 * phi_m.submit_value(scalar_product(avg_grad_pres, n_plus) - coef_jump*jump_pres, q);
1730 * phi_p.submit_gradient(-theta_p*0.5*jump_pres*n_plus, q);
1731 * phi_m.submit_gradient(-theta_p*0.5*jump_pres*n_plus, q);
1741 * The following function assembles boundary term
for the pressure
1747 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1748 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1752 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
1755 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
1758 *
if(boundary_id == 1) {
1760 * phi.gather_evaluate(src,
true,
true);
1762 *
const auto coef_jump = C_p*
std::abs((phi.get_normal_vector(0)*phi.inverse_jacobian(0))[dim - 1]);
1764 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1765 *
const auto& n_plus = phi.get_normal_vector(q);
1767 *
const auto& grad_pres = phi.get_gradient(q);
1768 *
const auto& pres = phi.get_value(q);
1770 * phi.submit_value(-scalar_product(grad_pres, n_plus) + coef_jump*pres , q);
1771 * phi.submit_normal_derivative(-theta_p*pres, q);
1781 * Before coding the
'apply_add' function, which is the one that will perform the
loop, we focus on
1782 * the linear system that arises to
project the
gradient of the pressure into the velocity space.
1783 * The following function assembles rhs cell term
for the projection of
gradient of pressure. Since no
1784 * integration by parts is performed, only a cell term contribution is present.
1790 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1791 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1795 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1801 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1802 * phi_pres.reinit(cell);
1807 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q)
1808 * phi.submit_value(phi_pres.get_gradient(q), q);
1817 * Put together all the previous steps
for porjection of pressure
gradient. Here we
loop only over cells
1823 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1824 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1825 * vmult_grad_p_projection(Vec& dst,
const Vec& src)
const {
1826 * this->data->
cell_loop(&NavierStokesProjectionOperator::assemble_rhs_cell_term_projection_grad_p,
1827 *
this, dst, src,
true);
1833 * Assemble now cell term
for the projection of
gradient of pressure. This is
nothing but a mass
matrix
1839 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1840 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1844 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1848 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1853 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q)
1854 * phi.submit_value(phi.get_value(q), q);
1863 * Put together all previous steps. This is the overriden function that effectively performs the
1864 *
matrix-vector multiplication.
1870 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1871 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1872 * apply_add(Vec& dst,
const Vec& src)
const {
1873 *
if(NS_stage == 1) {
1874 * this->data->
loop(&NavierStokesProjectionOperator::assemble_cell_term_velocity,
1875 * &NavierStokesProjectionOperator::assemble_face_term_velocity,
1876 * &NavierStokesProjectionOperator::assemble_boundary_term_velocity,
1877 *
this, dst, src,
false,
1881 *
else if(NS_stage == 2) {
1882 * this->data->
loop(&NavierStokesProjectionOperator::assemble_cell_term_pressure,
1883 * &NavierStokesProjectionOperator::assemble_face_term_pressure,
1884 * &NavierStokesProjectionOperator::assemble_boundary_term_pressure,
1885 *
this, dst, src,
false,
1889 *
else if(NS_stage == 3) {
1890 * this->data->
cell_loop(&NavierStokesProjectionOperator::assemble_cell_term_projection_grad_p,
1891 *
this, dst, src,
false);
1894 *
Assert(
false, ExcNotImplemented());
1900 * Finally, we focus on computing the
diagonal for preconditioners and we start by assembling
1901 * the
diagonal cell term
for the velocity. Since we
do not have access to the entries of the
matrix,
1902 * in order to compute the element i, we test the
matrix against a vector which is
equal to 1 in position i and 0 elsewhere.
1903 * This is why
'src' will result as unused.
1909 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
1910 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
1913 *
const unsigned int& ,
1914 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
1915 *
if(TR_BDF2_stage == 1) {
1917 * phi_old_extr(data, 0);
1923 *
for(
unsigned int d = 0;
d < dim; ++
d)
1924 * tmp[d] = make_vectorized_array<Number>(1.0);
1927 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1928 * phi_old_extr.reinit(cell);
1929 * phi_old_extr.gather_evaluate(u_extr,
true,
false);
1933 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
1934 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
1936 * phi.submit_dof_value(tmp, i);
1937 * phi.evaluate(
true,
true);
1940 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1941 *
const auto& u_int = phi.get_value(q);
1942 *
const auto& grad_u_int = phi.get_gradient(q);
1943 *
const auto& u_n_gamma_ov_2 = phi_old_extr.get_value(q);
1944 *
const auto& tensor_product_u_int =
outer_product(u_int, u_n_gamma_ov_2);
1946 * phi.submit_value(1.0/(gamma*dt)*u_int, q);
1947 * phi.submit_gradient(-a22*tensor_product_u_int + a22/Re*grad_u_int, q);
1949 * phi.integrate(
true,
true);
1950 *
diagonal[i] = phi.get_dof_value(i);
1952 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
1953 * phi.submit_dof_value(diagonal[i], i);
1954 * phi.distribute_local_to_global(dst);
1959 * phi_int_extr(data, 0);
1963 *
for(
unsigned int d = 0;
d < dim; ++
d)
1964 * tmp[d] = make_vectorized_array<Number>(1.0);
1967 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
1968 * phi_int_extr.reinit(cell);
1969 * phi_int_extr.gather_evaluate(u_extr,
true,
false);
1973 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
1974 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
1976 * phi.submit_dof_value(tmp, i);
1977 * phi.evaluate(
true,
true);
1980 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
1981 *
const auto& u_curr = phi.get_value(q);
1982 *
const auto& grad_u_curr = phi.get_gradient(q);
1983 *
const auto& u_n1_int = phi_int_extr.get_value(q);
1984 *
const auto& tensor_product_u_curr =
outer_product(u_curr, u_n1_int);
1986 * phi.submit_value(1.0/((1.0 - gamma)*dt)*u_curr, q);
1987 * phi.submit_gradient(-a33*tensor_product_u_curr + a33/Re*grad_u_curr, q);
1989 * phi.integrate(
true,
true);
1990 *
diagonal[i] = phi.get_dof_value(i);
1992 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
1993 * phi.submit_dof_value(diagonal[i], i);
1994 * phi.distribute_local_to_global(dst);
2002 * The following function assembles
diagonal face term
for the velocity
2008 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2009 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2012 *
const unsigned int& ,
2013 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
2014 *
if(TR_BDF2_stage == 1) {
2016 * phi_m(data,
false, 0),
2017 * phi_old_extr_p(data,
true, 0),
2018 * phi_old_extr_m(data,
false, 0);
2020 *
AssertDimension(phi_p.dofs_per_component, phi_m.dofs_per_component);
2024 * diagonal_m(phi_m.dofs_per_component);
2027 *
for(
unsigned int d = 0;
d < dim; ++
d)
2028 * tmp[d] = make_vectorized_array<Number>(1.0);
2031 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2032 * phi_old_extr_p.reinit(face);
2033 * phi_old_extr_p.gather_evaluate(u_extr,
true,
false);
2034 * phi_old_extr_m.reinit(face);
2035 * phi_old_extr_m.gather_evaluate(u_extr,
true,
false);
2036 * phi_p.reinit(face);
2037 * phi_m.reinit(face);
2039 *
const auto coef_jump = C_u*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
2040 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
2043 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2044 *
for(
unsigned int j = 0; j < phi_p.dofs_per_component; ++j) {
2048 * phi_p.submit_dof_value(tmp, i);
2049 * phi_p.evaluate(
true,
true);
2050 * phi_m.submit_dof_value(tmp, i);
2051 * phi_m.evaluate(
true,
true);
2054 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
2055 *
const auto& n_plus = phi_p.get_normal_vector(q);
2056 *
const auto& avg_grad_u_int = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
2057 *
const auto& jump_u_int = phi_p.get_value(q) - phi_m.get_value(q);
2058 *
const auto& avg_tensor_product_u_int = 0.5*(
outer_product(phi_p.get_value(q), phi_old_extr_p.get_value(q)) +
2059 *
outer_product(phi_m.get_value(q), phi_old_extr_m.get_value(q)));
2061 *
std::abs(scalar_product(phi_old_extr_m.get_value(q), n_plus)));
2063 * phi_p.submit_value(a22/Re*(-avg_grad_u_int*n_plus + coef_jump*jump_u_int) +
2064 * a22*avg_tensor_product_u_int*n_plus + 0.5*a22*lambda*jump_u_int , q);
2065 * phi_m.submit_value(-a22/Re*(-avg_grad_u_int*n_plus + coef_jump*jump_u_int) -
2066 * a22*avg_tensor_product_u_int*n_plus - 0.5*a22*lambda*jump_u_int, q);
2067 * phi_p.submit_normal_derivative(-theta_v*0.5*a22/Re*jump_u_int, q);
2068 * phi_m.submit_normal_derivative(-theta_v*0.5*a22/Re*jump_u_int, q);
2070 * phi_p.integrate(
true,
true);
2071 * diagonal_p[i] = phi_p.get_dof_value(i);
2072 * phi_m.integrate(
true,
true);
2073 * diagonal_m[i] = phi_m.get_dof_value(i);
2075 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2076 * phi_p.submit_dof_value(diagonal_p[i], i);
2077 * phi_m.submit_dof_value(diagonal_m[i], i);
2079 * phi_p.distribute_local_to_global(dst);
2080 * phi_m.distribute_local_to_global(dst);
2085 * phi_m(data,
false, 0),
2086 * phi_extr_p(data,
true, 0),
2087 * phi_extr_m(data,
false, 0);
2089 *
AssertDimension(phi_p.dofs_per_component, phi_m.dofs_per_component);
2091 * diagonal_m(phi_m.dofs_per_component);
2093 *
for(
unsigned int d = 0;
d < dim; ++
d)
2094 * tmp[d] = make_vectorized_array<Number>(1.0);
2097 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2098 * phi_extr_p.reinit(face);
2099 * phi_extr_p.gather_evaluate(u_extr,
true,
false);
2100 * phi_extr_m.reinit(face);
2101 * phi_extr_m.gather_evaluate(u_extr,
true,
false);
2102 * phi_p.reinit(face);
2103 * phi_m.reinit(face);
2105 *
const auto coef_jump = C_u*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
2106 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
2109 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2110 *
for(
unsigned int j = 0; j < phi_p.dofs_per_component; ++j) {
2114 * phi_p.submit_dof_value(tmp, i);
2115 * phi_p.evaluate(
true,
true);
2116 * phi_m.submit_dof_value(tmp, i);
2117 * phi_m.evaluate(
true,
true);
2120 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
2121 *
const auto& n_plus = phi_p.get_normal_vector(q);
2122 *
const auto& avg_grad_u = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
2123 *
const auto& jump_u = phi_p.get_value(q) - phi_m.get_value(q);
2124 *
const auto& avg_tensor_product_u = 0.5*(
outer_product(phi_p.get_value(q), phi_extr_p.get_value(q)) +
2125 *
outer_product(phi_m.get_value(q), phi_extr_m.get_value(q)));
2127 *
std::abs(scalar_product(phi_extr_m.get_value(q), n_plus)));
2129 * phi_p.submit_value(a33/Re*(-avg_grad_u*n_plus + coef_jump*jump_u) +
2130 * a33*avg_tensor_product_u*n_plus + 0.5*a33*lambda*jump_u, q);
2131 * phi_m.submit_value(-a33/Re*(-avg_grad_u*n_plus + coef_jump*jump_u) -
2132 * a33*avg_tensor_product_u*n_plus - 0.5*a33*lambda*jump_u, q);
2133 * phi_p.submit_normal_derivative(-theta_v*0.5*a33/Re*jump_u, q);
2134 * phi_m.submit_normal_derivative(-theta_v*0.5*a33/Re*jump_u, q);
2136 * phi_p.integrate(
true,
true);
2137 * diagonal_p[i] = phi_p.get_dof_value(i);
2138 * phi_m.integrate(
true,
true);
2139 * diagonal_m[i] = phi_m.get_dof_value(i);
2141 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2142 * phi_p.submit_dof_value(diagonal_p[i], i);
2143 * phi_m.submit_dof_value(diagonal_m[i], i);
2145 * phi_p.distribute_local_to_global(dst);
2146 * phi_m.distribute_local_to_global(dst);
2154 * The following function assembles boundary term
for the velocity
2160 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2161 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2164 *
const unsigned int& ,
2165 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
2166 *
if(TR_BDF2_stage == 1) {
2168 * phi_old_extr(data,
true, 0);
2172 *
for(
unsigned int d = 0;
d < dim; ++
d)
2173 * tmp[d] = make_vectorized_array<Number>(1.0);
2176 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2177 * phi_old_extr.reinit(face);
2178 * phi_old_extr.gather_evaluate(u_extr,
true,
false);
2182 *
const auto coef_jump = C_u*
std::abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
2184 *
if(boundary_id != 1) {
2185 *
const double coef_trasp = 0.0;
2188 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2189 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2191 * phi.submit_dof_value(tmp, i);
2192 * phi.evaluate(
true,
true);
2195 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2196 *
const auto& n_plus = phi.get_normal_vector(q);
2197 *
const auto& grad_u_int = phi.get_gradient(q);
2198 *
const auto& u_int = phi.get_value(q);
2199 *
const auto& tensor_product_u_int =
outer_product(phi.get_value(q), phi_old_extr.get_value(q));
2200 *
const auto&
lambda =
std::abs(scalar_product(phi_old_extr.get_value(q), n_plus));
2202 * phi.submit_value(a22/Re*(-grad_u_int*n_plus + 2.0*coef_jump*u_int) +
2203 * a22*coef_trasp*tensor_product_u_int*n_plus + a22*lambda*u_int, q);
2204 * phi.submit_normal_derivative(-theta_v*a22/Re*u_int, q);
2206 * phi.integrate(
true,
true);
2207 *
diagonal[i] = phi.get_dof_value(i);
2209 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2210 * phi.submit_dof_value(diagonal[i], i);
2211 * phi.distribute_local_to_global(dst);
2215 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2216 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2218 * phi.submit_dof_value(tmp, i);
2219 * phi.evaluate(
true,
true);
2222 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2223 *
const auto& n_plus = phi.get_normal_vector(q);
2224 *
const auto& grad_u_int = phi.get_gradient(q);
2225 *
const auto& u_int = phi.get_value(q);
2226 *
const auto&
lambda =
std::abs(scalar_product(phi_old_extr.get_value(q), n_plus));
2228 *
const auto& point_vectorized = phi.quadrature_point(q);
2229 *
auto u_int_m = u_int;
2230 *
auto grad_u_int_m = grad_u_int;
2231 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
2233 *
for(
unsigned int d = 0;
d < dim; ++
d)
2234 * point[d] = point_vectorized[d][v];
2236 * u_int_m[1][v] = -u_int_m[1][v];
2238 * grad_u_int_m[0][0][v] = -grad_u_int_m[0][0][v];
2239 * grad_u_int_m[0][1][v] = -grad_u_int_m[0][1][v];
2242 * phi.submit_value(a22/Re*(-(0.5*(grad_u_int + grad_u_int_m))*n_plus + coef_jump*(u_int - u_int_m)) +
2243 * a22*
outer_product(0.5*(u_int + u_int_m), phi_old_extr.get_value(q))*n_plus +
2244 * a22*0.5*
lambda*(u_int - u_int_m), q);
2245 * phi.submit_normal_derivative(-theta_v*a22/Re*(u_int - u_int_m), q);
2247 * phi.integrate(
true,
true);
2248 *
diagonal[i] = phi.get_dof_value(i);
2250 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2251 * phi.submit_dof_value(diagonal[i], i);
2252 * phi.distribute_local_to_global(dst);
2258 * phi_extr(data,
true, 0);
2262 *
for(
unsigned int d = 0;
d < dim; ++
d)
2263 * tmp[d] = make_vectorized_array<Number>(1.0);
2266 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2267 * phi_extr.reinit(face);
2268 * phi_extr.gather_evaluate(u_extr,
true,
false);
2272 *
const auto coef_jump = C_u*
std::abs((phi.get_normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]);
2274 *
if(boundary_id != 1) {
2275 *
const double coef_trasp = 0.0;
2278 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2279 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2281 * phi.submit_dof_value(tmp, i);
2282 * phi.evaluate(
true,
true);
2285 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2286 *
const auto& n_plus = phi.get_normal_vector(q);
2287 *
const auto& grad_u = phi.get_gradient(q);
2288 *
const auto& u = phi.get_value(q);
2289 *
const auto& tensor_product_u =
outer_product(phi.get_value(q), phi_extr.get_value(q));
2290 *
const auto&
lambda =
std::abs(scalar_product(phi_extr.get_value(q), n_plus));
2292 * phi.submit_value(a33/Re*(-grad_u*n_plus + 2.0*coef_jump*u) +
2293 * a33*coef_trasp*tensor_product_u*n_plus + a33*lambda*u, q);
2294 * phi.submit_normal_derivative(-theta_v*a33/Re*u, q);
2296 * phi.integrate(
true,
true);
2297 *
diagonal[i] = phi.get_dof_value(i);
2299 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2300 * phi.submit_dof_value(diagonal[i], i);
2301 * phi.distribute_local_to_global(dst);
2305 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2306 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2308 * phi.submit_dof_value(tmp, i);
2309 * phi.evaluate(
true,
true);
2312 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2313 *
const auto& n_plus = phi.get_normal_vector(q);
2314 *
const auto& grad_u = phi.get_gradient(q);
2315 *
const auto& u = phi.get_value(q);
2316 *
const auto&
lambda =
std::abs(scalar_product(phi_extr.get_value(q), n_plus));
2318 *
const auto& point_vectorized = phi.quadrature_point(q);
2320 *
auto grad_u_m = grad_u;
2321 *
for(
unsigned int v = 0; v < VectorizedArray<Number>::size(); ++v) {
2323 *
for(
unsigned int d = 0;
d < dim; ++
d)
2324 * point[d] = point_vectorized[d][v];
2326 * u_m[1][v] = -u_m[1][v];
2328 * grad_u_m[0][0][v] = -grad_u_m[0][0][v];
2329 * grad_u_m[0][1][v] = -grad_u_m[0][1][v];
2332 * phi.submit_value(a33/Re*(-(0.5*(grad_u + grad_u_m))*n_plus + coef_jump*(u - u_m)) +
2333 * a33*
outer_product(0.5*(u + u_m), phi_extr.get_value(q))*n_plus +
2334 * a33*0.5*
lambda*(u - u_m), q);
2335 * phi.submit_normal_derivative(-theta_v*a33/Re*(u - u_m), q);
2337 * phi.integrate(
true,
true);
2338 *
diagonal[i] = phi.get_dof_value(i);
2340 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2341 * phi.submit_dof_value(diagonal[i], i);
2342 * phi.distribute_local_to_global(dst);
2351 * Now we consider the pressure related bilinear forms. We
first assemble diagonal cell term
for the pressure
2357 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2358 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2361 *
const unsigned int& ,
2362 *
const std::pair<unsigned int, unsigned int>& cell_range)
const {
2369 *
const double coeff = (TR_BDF2_stage == 1) ? 1e6*gamma*dt*gamma*dt : 1e6*(1.0 -
gamma)*dt*(1.0 -
gamma)*dt;
2372 *
for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell) {
2376 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2377 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2379 * phi.submit_dof_value(make_vectorized_array<Number>(1.0), i);
2385 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2386 * phi.submit_value(1.0/coeff*phi.get_value(q), q);
2387 * phi.submit_gradient(phi.get_gradient(q), q);
2390 *
diagonal[i] = phi.get_dof_value(i);
2392 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2393 * phi.submit_dof_value(diagonal[i], i);
2395 * phi.distribute_local_to_global(dst);
2402 * The following function assembles
diagonal face term
for the pressure
2408 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2409 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2412 *
const unsigned int& ,
2413 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
2415 * phi_m(data,
false, 1, 1);
2417 *
AssertDimension(phi_p.dofs_per_component, phi_m.dofs_per_component);
2419 * diagonal_m(phi_m.dofs_per_component);
2424 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2425 * phi_p.reinit(face);
2426 * phi_m.reinit(face);
2428 *
const auto coef_jump = C_p*0.5*(
std::abs((phi_p.get_normal_vector(0)*phi_p.inverse_jacobian(0))[dim - 1]) +
2429 *
std::abs((phi_m.get_normal_vector(0)*phi_m.inverse_jacobian(0))[dim - 1]));
2432 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2433 *
for(
unsigned int j = 0; j < phi_p.dofs_per_component; ++j) {
2437 * phi_p.submit_dof_value(make_vectorized_array<Number>(1.0), i);
2438 * phi_m.submit_dof_value(make_vectorized_array<Number>(1.0), i);
2443 *
for(
unsigned int q = 0; q < phi_p.n_q_points; ++q) {
2444 *
const auto& n_plus = phi_p.get_normal_vector(q);
2446 *
const auto& avg_grad_pres = 0.5*(phi_p.get_gradient(q) + phi_m.get_gradient(q));
2447 *
const auto& jump_pres = phi_p.get_value(q) - phi_m.get_value(q);
2449 * phi_p.submit_value(-scalar_product(avg_grad_pres, n_plus) + coef_jump*jump_pres, q);
2450 * phi_m.submit_value(scalar_product(avg_grad_pres, n_plus) - coef_jump*jump_pres, q);
2451 * phi_p.submit_gradient(-theta_p*0.5*jump_pres*n_plus, q);
2452 * phi_m.submit_gradient(-theta_p*0.5*jump_pres*n_plus, q);
2455 * diagonal_p[i] = phi_p.get_dof_value(i);
2457 * diagonal_m[i] = phi_m.get_dof_value(i);
2459 *
for(
unsigned int i = 0; i < phi_p.dofs_per_component; ++i) {
2460 * phi_p.submit_dof_value(diagonal_p[i], i);
2461 * phi_m.submit_dof_value(diagonal_m[i], i);
2463 * phi_p.distribute_local_to_global(dst);
2464 * phi_m.distribute_local_to_global(dst);
2477 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2478 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2481 *
const unsigned int& ,
2482 *
const std::pair<unsigned int, unsigned int>& face_range)
const {
2487 *
for(
unsigned int face = face_range.first; face < face_range.second; ++face) {
2490 *
if(boundary_id == 1) {
2493 *
const auto coef_jump = C_p*
std::abs((phi.get_normal_vector(0)*phi.inverse_jacobian(0))[dim - 1]);
2495 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i) {
2496 *
for(
unsigned int j = 0; j < phi.dofs_per_component; ++j)
2498 * phi.submit_dof_value(make_vectorized_array<Number>(1.0), i);
2501 *
for(
unsigned int q = 0; q < phi.n_q_points; ++q) {
2502 *
const auto& n_plus = phi.get_normal_vector(q);
2504 *
const auto& grad_pres = phi.get_gradient(q);
2505 *
const auto& pres = phi.get_value(q);
2507 * phi.submit_value(-scalar_product(grad_pres, n_plus) + 2.0*coef_jump*pres , q);
2508 * phi.submit_normal_derivative(-theta_p*pres, q);
2511 *
diagonal[i] = phi.get_dof_value(i);
2513 *
for(
unsigned int i = 0; i < phi.dofs_per_component; ++i)
2514 * phi.submit_dof_value(diagonal[i], i);
2515 * phi.distribute_local_to_global(dst);
2523 * Put together all previous steps. We create a dummy auxliary vector that serves
for the src input argument in
2524 * the previous
functions that as we have seen before is unused. Then everything is done by the
'loop' function
2525 * and it is saved in the field
'inverse_diagonal_entries' already present in the base
class. Anyway since there is
2526 * only one field, we need to resize properly depending on whether we are considering the velocity or the pressure.
2532 *
template<
int dim,
int fe_degree_p,
int fe_degree_v,
int n_q_po
ints_1d_p,
int n_q_po
ints_1d_v,
typename Vec>
2533 *
void NavierStokesProjectionOperator<dim, fe_degree_p, fe_degree_v, n_q_points_1d_p, n_q_points_1d_v, Vec>::
2535 *
Assert(NS_stage == 1 || NS_stage == 2, ExcInternalError());
2538 *
auto& inverse_diagonal = this->inverse_diagonal_entries->get_vector();
2540 *
if(NS_stage == 1) {
2541 * ::MatrixFreeTools::compute_diagonal<dim, Number, VectorizedArray<Number>>
2544 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& cell_range) {
2545 * (this->assemble_diagonal_cell_term_velocity)(data, dst, src, cell_range);
2547 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& face_range) {
2548 * (this->assemble_diagonal_face_term_velocity)(data, dst, src, face_range);
2550 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& boundary_range) {
2551 * (this->assemble_diagonal_boundary_term_velocity)(data, dst, src, boundary_range);
2555 *
else if(NS_stage == 2) {
2556 * ::MatrixFreeTools::compute_diagonal<dim, Number, VectorizedArray<Number>>
2559 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& cell_range) {
2560 * (this->assemble_diagonal_cell_term_pressure)(data, dst, src, cell_range);
2562 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& face_range) {
2563 * (this->assemble_diagonal_face_term_pressure)(data, dst, src, face_range);
2565 * [&](
const auto& data,
auto& dst,
const auto& src,
const auto& boundary_range) {
2566 * (this->assemble_diagonal_boundary_term_pressure)(data, dst, src, boundary_range);
2571 *
for(
unsigned int i = 0; i < inverse_diagonal.locally_owned_size(); ++i) {
2572 *
Assert(inverse_diagonal.local_element(i) != 0.0,
2573 * ExcMessage(
"No diagonal entry in a definite operator should be zero"));
2574 * inverse_diagonal.local_element(i) = 1.0/inverse_diagonal.local_element(i);
2583 * @sect{The <code>NavierStokesProjection</code>
class}
2587 * Now we are ready
for the main
class of the program. It
implements the calls to the various steps
2588 * of the projection method for Navier-Stokes equations.
2595 * class NavierStokesProjection {
2597 * NavierStokesProjection(RunTimeParameters::Data_Storage& data);
2599 *
void run(
const bool verbose =
false,
const unsigned int output_interval = 10);
2604 *
const double gamma;
2605 *
unsigned int TR_BDF2_stage;
2609 * EquationData::Velocity<dim> vel_init;
2610 * EquationData::Pressure<dim> pres_init;
2648 * <<
" The time step " << arg1 <<
" is out of range."
2650 * <<
" The permitted range is (0," << arg2 <<
"]");
2654 *
void setup_dofs();
2656 *
void initialize();
2658 *
void interpolate_velocity();
2660 *
void diffusion_step();
2662 *
void projection_step();
2664 *
void project_grad(
const unsigned int flag);
2666 *
double get_maximal_velocity();
2668 *
double get_maximal_difference_velocity();
2670 *
void output_results(
const unsigned int step);
2672 *
void refine_mesh();
2674 *
void interpolate_max_res(
const unsigned int level);
2676 *
void save_max_res();
2679 *
void compute_lift_and_drag();
2682 * std::shared_ptr<MatrixFree<dim, double>> matrix_free_storage;
2685 * NavierStokesProjectionOperator<dim, EquationData::degree_p, EquationData::degree_p + 1,
2686 * EquationData::degree_p + 1, EquationData::degree_p + 2,
2690 *
MGLevelObject<NavierStokesProjectionOperator<dim, EquationData::degree_p, EquationData::degree_p + 1,
2691 * EquationData::degree_p + 1, EquationData::degree_p + 2,
2700 * constraints_pressure;
2703 *
unsigned int max_its;
2706 *
unsigned int max_loc_refinements;
2707 *
unsigned int min_loc_refinements;
2708 *
unsigned int refinement_iterations;
2710 * std::string saving_dir;
2715 * std::ofstream time_out;
2719 * std::ofstream output_n_dofs_velocity;
2720 * std::ofstream output_n_dofs_pressure;
2722 * std::ofstream output_lift;
2723 * std::ofstream output_drag;
2729 * In the constructor, we just read all the data from the
2730 * <code>Data_Storage</code>
object that is passed as an argument, verify that
2731 * the data we read are reasonable and,
finally, create the
triangulation and
2739 * NavierStokesProjection<dim>::NavierStokesProjection(RunTimeParameters::Data_Storage& data):
2740 * t_0(data.initial_time),
2741 * T(data.final_time),
2744 * Re(data.Reynolds),
2746 * vel_init(data.initial_time),
2747 * pres_init(data.initial_time),
2750 * fe_velocity(
FE_DGQ<dim>(EquationData::degree_p + 1), dim),
2751 * fe_pressure(
FE_DGQ<dim>(EquationData::degree_p), 1),
2754 * quadrature_pressure(EquationData::degree_p + 1),
2755 * quadrature_velocity(EquationData::degree_p + 2),
2756 * navier_stokes_matrix(data),
2757 * max_its(data.max_iterations),
2759 * max_loc_refinements(data.max_loc_refinements),
2760 * min_loc_refinements(data.min_loc_refinements),
2761 * refinement_iterations(data.refinement_iterations),
2762 * saving_dir(data.dir),
2764 * time_out(
"./" + data.dir +
"/time_analysis_" +
2768 * output_n_dofs_velocity(
"./" + data.dir +
"/n_dofs_velocity.dat",
std::ofstream::out),
2769 * output_n_dofs_pressure(
"./" + data.dir +
"/n_dofs_pressure.dat",
std::ofstream::out),
2770 * output_lift(
"./" + data.dir +
"/lift.dat",
std::ofstream::out),
2771 * output_drag(
"./" + data.dir +
"/drag.dat",
std::ofstream::out) {
2772 *
if(EquationData::degree_p < 1) {
2774 * <<
" WARNING: The chosen pair of finite element spaces is not stable."
2776 * <<
" The obtained results will be nonsense" << std::endl;
2779 *
AssertThrow(!((dt <= 0.0) || (dt > 0.5*T)), ExcInvalidTimeStep(dt, 0.5*T));
2781 * matrix_free_storage = std::make_shared<MatrixFree<dim, double>>();
2791 * The method that creates the
triangulation and refines it the needed number
2799 *
void NavierStokesProjection<dim>::create_triangulation(
const unsigned int n_refines) {
2802 *
GridGenerator::plate_with_a_hole(
triangulation, 0.5, 1.0, 1.0, 1.1, 1.0, 19.0,
Point<2>(2.0, 2.0), 0, 1, 1.0, 2,
true);
2805 * pcout <<
"Number of refines = " << n_refines << std::endl;
2812 * After creating the
triangulation, it creates the mesh dependent
2813 * data, i.e. it distributes degrees of freedom, and
2814 * initializes the vectors that we will use.
2821 *
void NavierStokesProjection<dim>::setup_dofs() {
2822 * pcout <<
"Number of active cells: " <<
triangulation.n_global_active_cells() << std::endl;
2823 * pcout <<
"Number of levels: " <<
triangulation.n_global_levels() << std::endl;
2826 * dof_handler_velocity.distribute_dofs(fe_velocity);
2827 * dof_handler_pressure.distribute_dofs(fe_pressure);
2829 * pcout <<
"dim (X_h) = " << dof_handler_velocity.n_dofs()
2831 * <<
"dim (M_h) = " << dof_handler_pressure.n_dofs()
2833 * <<
"Re = " << Re << std::endl
2837 * output_n_dofs_velocity << dof_handler_velocity.n_dofs() << std::endl;
2838 * output_n_dofs_pressure << dof_handler_pressure.n_dofs() << std::endl;
2850 * std::vector<const DoFHandler<dim>*> dof_handlers;
2853 * dof_handlers.push_back(&dof_handler_velocity);
2854 * dof_handlers.push_back(&dof_handler_pressure);
2856 * constraints_velocity.
clear();
2857 * constraints_velocity.close();
2858 * constraints_pressure.clear();
2859 * constraints_pressure.close();
2860 * std::vector<const AffineConstraints<double>*> constraints;
2861 * constraints.push_back(&constraints_velocity);
2862 * constraints.push_back(&constraints_pressure);
2864 * std::vector<QGauss<1>> quadratures;
2868 * quadratures.push_back(
QGauss<1>(EquationData::degree_p + 2));
2869 * quadratures.push_back(
QGauss<1>(EquationData::degree_p + 1));
2873 * matrix_free_storage->reinit(
MappingQ1<dim>(),dof_handlers, constraints, quadratures, additional_data);
2874 * matrix_free_storage->initialize_dof_vector(u_star, 0);
2875 * matrix_free_storage->initialize_dof_vector(rhs_u, 0);
2876 * matrix_free_storage->initialize_dof_vector(u_n, 0);
2877 * matrix_free_storage->initialize_dof_vector(u_extr, 0);
2878 * matrix_free_storage->initialize_dof_vector(u_n_minus_1, 0);
2879 * matrix_free_storage->initialize_dof_vector(u_n_gamma, 0);
2880 * matrix_free_storage->initialize_dof_vector(u_tmp, 0);
2881 * matrix_free_storage->initialize_dof_vector(grad_pres_int, 0);
2883 * matrix_free_storage->initialize_dof_vector(pres_int, 1);
2884 * matrix_free_storage->initialize_dof_vector(pres_n, 1);
2885 * matrix_free_storage->initialize_dof_vector(rhs_p, 1);
2891 * mg_matrices.clear_elements();
2892 * dof_handler_velocity.distribute_mg_dofs();
2893 * dof_handler_pressure.distribute_mg_dofs();
2895 *
const unsigned int nlevels =
triangulation.n_global_levels();
2896 * mg_matrices.resize(0, nlevels - 1);
2905 * std::vector<const DoFHandler<dim>*> dof_handlers_mg;
2906 * dof_handlers_mg.push_back(&dof_handler_velocity);
2907 * dof_handlers_mg.push_back(&dof_handler_pressure);
2908 * std::vector<const AffineConstraints<float>*> constraints_mg;
2910 * constraints_velocity_mg.
clear();
2911 * constraints_velocity_mg.close();
2912 * constraints_mg.push_back(&constraints_velocity_mg);
2914 * constraints_pressure_mg.
clear();
2915 * constraints_pressure_mg.close();
2916 * constraints_mg.push_back(&constraints_pressure_mg);
2919 * mg_mf_storage_level->reinit(
MappingQ1<dim>(),dof_handlers_mg, constraints_mg, quadratures, additional_data_mg);
2920 *
const std::vector<unsigned int> tmp = {1};
2921 * mg_matrices[
level].initialize(mg_mf_storage_level, tmp, tmp);
2922 * mg_matrices[
level].set_dt(dt);
2923 * mg_matrices[
level].set_NS_stage(2);
2926 * Linfty_error_per_cell_vel.reinit(
triangulation.n_active_cells());
2932 * This method loads the
initial data. It simply uses the class <code>Pressure</code> instance
for the pressure
2933 * and the class <code>Velocity</code> instance
for the velocity.
2940 *
void NavierStokesProjection<dim>::initialize() {
2952 * This function computes the extrapolated velocity to be used in the momentum predictor
2959 *
void NavierStokesProjection<dim>::interpolate_velocity() {
2964 * --- TR-BDF2
first step
2967 *
if(TR_BDF2_stage == 1) {
2968 * u_extr.equ(1.0 + gamma/(2.0*(1.0 - gamma)), u_n);
2969 * u_tmp.equ(gamma/(2.0*(1.0 - gamma)), u_n_minus_1);
2974 * --- TR-BDF2
second step
2978 * u_extr.equ(1.0 + (1.0 - gamma)/gamma, u_n_gamma);
2979 * u_tmp.equ((1.0 - gamma)/gamma, u_n);
2987 * We are
finally ready to solve the diffusion step.
2994 *
void NavierStokesProjection<dim>::diffusion_step() {
2999 *
const std::vector<unsigned int> tmp = {0};
3000 * navier_stokes_matrix.initialize(matrix_free_storage, tmp, tmp);
3003 * navier_stokes_matrix.set_NS_stage(1);
3008 *
if(TR_BDF2_stage == 1) {
3009 * navier_stokes_matrix.vmult_rhs_velocity(rhs_u, {u_n, u_extr, pres_n});
3010 * navier_stokes_matrix.set_u_extr(u_extr);
3014 * navier_stokes_matrix.vmult_rhs_velocity(rhs_u, {u_n, u_n_gamma, pres_int, u_extr});
3015 * navier_stokes_matrix.set_u_extr(u_extr);
3020 *
SolverControl solver_control(max_its, eps*rhs_u.l2_norm());
3025 * EquationData::degree_p,
3026 * EquationData::degree_p + 1,
3027 * EquationData::degree_p + 1,
3028 * EquationData::degree_p + 2,
3030 * navier_stokes_matrix.compute_diagonal();
3031 * preconditioner.initialize(navier_stokes_matrix);
3033 * gmres.solve(navier_stokes_matrix, u_star, rhs_u, preconditioner);
3039 * Next, we solve the projection step.
3046 *
void NavierStokesProjection<dim>::projection_step() {
3051 *
const std::vector<unsigned int> tmp = {1};
3052 * navier_stokes_matrix.initialize(matrix_free_storage, tmp, tmp);
3054 * navier_stokes_matrix.set_NS_stage(2);
3056 *
if(TR_BDF2_stage == 1)
3057 * navier_stokes_matrix.vmult_rhs_pressure(rhs_p, {u_star, pres_n});
3059 * navier_stokes_matrix.vmult_rhs_pressure(rhs_p, {u_star, pres_int});
3062 *
SolverControl solver_control(max_its, eps*rhs_p.l2_norm());
3067 * mg_transfer.
build(dof_handler_pressure);
3070 * EquationData::degree_p,
3071 * EquationData::degree_p + 1,
3072 * EquationData::degree_p + 1,
3073 * EquationData::degree_p + 2,
3081 * smoother_data[
level].smoothing_range = 15.0;
3082 * smoother_data[
level].degree = 3;
3083 * smoother_data[
level].eig_cg_n_iterations = 10;
3086 * smoother_data[0].smoothing_range = 2
e-2;
3088 * smoother_data[0].eig_cg_n_iterations = mg_matrices[0].m();
3090 * mg_matrices[
level].compute_diagonal();
3091 * smoother_data[
level].preconditioner = mg_matrices[
level].get_matrix_diagonal_inverse();
3093 * mg_smoother.initialize(mg_matrices, smoother_data);
3099 * NavierStokesProjectionOperator<dim,
3100 * EquationData::degree_p,
3101 * EquationData::degree_p + 1,
3102 * EquationData::degree_p + 1,
3103 * EquationData::degree_p + 2,
3116 *
if(TR_BDF2_stage == 1) {
3117 * pres_int = pres_n;
3118 * cg.solve(navier_stokes_matrix, pres_int, rhs_p, preconditioner);
3121 * pres_n = pres_int;
3122 * cg.solve(navier_stokes_matrix, pres_n, rhs_p, preconditioner);
3129 * This implements the projection step
for the
gradient of pressure
3136 *
void NavierStokesProjection<dim>::project_grad(
const unsigned int flag) {
3141 *
Assert(flag > 0, ExcInternalError());
3144 *
const std::vector<unsigned int> tmp = {0};
3145 * navier_stokes_matrix.initialize(matrix_free_storage, tmp, tmp);
3148 * navier_stokes_matrix.vmult_grad_p_projection(rhs_u, pres_n);
3149 *
else if(flag == 2)
3150 * navier_stokes_matrix.vmult_grad_p_projection(rhs_u, pres_int);
3153 * navier_stokes_matrix.set_NS_stage(3);
3156 *
SolverControl solver_control(max_its, 1e-12*rhs_u.l2_norm());
3164 * The following function is used in determining the maximal velocity
3165 * in order to compute the Courant number.
3172 *
double NavierStokesProjection<dim>::get_maximal_velocity() {
3173 *
return u_n.linfty_norm();
3179 * The following function is used in determining the maximal nodal difference
3180 * between old and current velocity
value in order to see
if we have reched steady-state.
3187 *
double NavierStokesProjection<dim>::get_maximal_difference_velocity() {
3189 * u_tmp -= u_n_minus_1;
3191 *
return u_tmp.linfty_norm();
3197 * This method plots the current solution. The main difficulty is that we want
3198 * to create a single output file that contains the data
for all velocity
3199 * components and the pressure. On the other hand, velocities and the pressure
3200 * live on separate
DoFHandler objects, so we need to pay attention when we use
3201 *
'add_data_vector' to select the proper space.
3208 *
void NavierStokesProjection<dim>::output_results(
const unsigned int step) {
3213 * std::vector<std::string> velocity_names(dim,
"v");
3214 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
3216 * u_n.update_ghost_values();
3217 * data_out.add_data_vector(dof_handler_velocity, u_n, velocity_names, component_interpretation_velocity);
3218 * pres_n.update_ghost_values();
3221 * std::vector<std::string> velocity_names_old(dim,
"v_old");
3222 * u_n_minus_1.update_ghost_values();
3223 * data_out.add_data_vector(dof_handler_velocity, u_n_minus_1, velocity_names_old, component_interpretation_velocity);
3226 * PostprocessorVorticity<dim> postprocessor;
3227 * data_out.add_data_vector(dof_handler_velocity, u_n, postprocessor);
3232 * data_out.write_vtu_in_parallel(output, MPI_COMM_WORLD);
3240 * @sect{<code>NavierStokesProjection::compute_lift_and_drag</code>}
3244 * This routine computes the lift and the drag forces in a non-dimensional framework
3245 * (so basically
for the classical coefficients, it is necessary to multiply by a factor 2).
3252 *
void NavierStokesProjection<dim>::compute_lift_and_drag() {
3253 *
QGauss<dim - 1> face_quadrature_formula(EquationData::degree_p + 2);
3254 *
const int n_q_points = face_quadrature_formula.size();
3256 * std::vector<double> pressure_values(n_q_points);
3257 * std::vector<std::vector<Tensor<1, dim>>> velocity_gradients(n_q_points, std::vector<
Tensor<1, dim>>(dim));
3266 *
FEFaceValues<dim> fe_face_values_velocity(fe_velocity, face_quadrature_formula,
3271 *
double local_drag = 0.0;
3272 *
double local_lift = 0.0;
3278 *
auto tmp_cell = dof_handler_pressure.begin_active();
3279 *
for(
const auto& cell : dof_handler_velocity.active_cell_iterators()) {
3280 *
if(cell->is_locally_owned()) {
3281 *
for(
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face) {
3282 *
if(cell->face(face)->at_boundary() && cell->face(face)->boundary_id() == 4) {
3283 * fe_face_values_velocity.reinit(cell, face);
3284 * fe_face_values_pressure.reinit(tmp_cell, face);
3286 * fe_face_values_velocity.get_function_gradients(u_n, velocity_gradients);
3287 * fe_face_values_pressure.get_function_values(pres_n, pressure_values);
3289 *
for(
int q = 0; q < n_q_points; q++) {
3290 * normal_vector = -fe_face_values_velocity.normal_vector(q);
3292 *
for(
unsigned int d = 0;
d < dim; ++
d) {
3293 * fluid_pressure[
d][
d] = pressure_values[q];
3294 *
for(
unsigned int k = 0; k < dim; ++k)
3295 * fluid_stress[d][k] = 1.0/Re*velocity_gradients[q][d][k];
3297 * fluid_stress = fluid_stress - fluid_pressure;
3299 * forces = fluid_stress*normal_vector*fe_face_values_velocity.JxW(q);
3301 * local_drag += forces[0];
3302 * local_lift += forces[1];
3315 * output_lift << lift << std::endl;
3316 * output_drag << drag << std::endl;
3325 * @sect{ <code>NavierStokesProjection::refine_mesh</code>}
3329 * After finding a good
initial guess on the coarse mesh, we hope to
3330 * decrease the error through refining the mesh. We also need to transfer the current solution to the
3337 *
template <
int dim>
3338 *
void NavierStokesProjection<dim>::refine_mesh() {
3345 * tmp_velocity.
reinit(dof_handler_velocity.locally_owned_dofs(), locally_relevant_dofs, MPI_COMM_WORLD);
3346 * tmp_velocity = u_n;
3347 * tmp_velocity.update_ghost_values();
3354 *
const auto cell_worker = [&](
const Iterator& cell,
3355 * ScratchData<dim>& scratch_data,
3356 * CopyData& copy_data) {
3358 * fe_values.
reinit(cell);
3361 * std::vector<std::vector<Tensor<1, dim>>>
gradients(fe_values.n_quadrature_points, std::vector<
Tensor<1, dim>>(dim));
3362 * fe_values.get_function_gradients(tmp_velocity, gradients);
3363 * copy_data.cell_index = cell->active_cell_index();
3364 *
double vorticity_norm_square = 0.0;
3367 *
for(
unsigned k = 0; k < fe_values.n_quadrature_points; ++k) {
3369 * vorticity_norm_square += vorticity*vorticity*fe_values.JxW(k);
3371 * copy_data.value = cell->diameter()*cell->diameter()*vorticity_norm_square;
3376 *
const auto copier = [&](
const CopyData ©_data) {
3378 * estimated_error_per_cell[copy_data.cell_index] += copy_data.value;
3382 * ScratchData<dim> scratch_data(fe_velocity, EquationData::degree_p + 2, cell_flags);
3383 * CopyData copy_data;
3385 * dof_handler_velocity.end(),
3395 *
for(
const auto& cell:
triangulation.active_cell_iterators()) {
3396 *
if(cell->refine_flag_set() &&
static_cast<unsigned int>(cell->level()) == max_loc_refinements)
3397 * cell->clear_refine_flag();
3398 *
if(cell->coarsen_flag_set() &&
static_cast<unsigned int>(cell->level()) == min_loc_refinements)
3399 * cell->clear_coarsen_flag();
3406 * std::vector<const LinearAlgebra::distributed::Vector<double>*> velocities;
3407 * velocities.push_back(&u_n);
3408 * velocities.push_back(&u_n_minus_1);
3410 * solution_transfer_velocity(dof_handler_velocity);
3411 * solution_transfer_velocity.prepare_for_coarsening_and_refinement(velocities);
3413 * solution_transfer_pressure(dof_handler_pressure);
3414 * solution_transfer_pressure.prepare_for_coarsening_and_refinement(pres_n);
3425 * transfer_velocity_minus_1,
3426 * transfer_pressure;
3427 * transfer_velocity.
reinit(u_n);
3428 * transfer_velocity.zero_out_ghost_values();
3429 * transfer_velocity_minus_1.reinit(u_n_minus_1);
3430 * transfer_velocity_minus_1.zero_out_ghost_values();
3431 * transfer_pressure.reinit(pres_n);
3432 * transfer_pressure.zero_out_ghost_values();
3434 * std::vector<LinearAlgebra::distributed::Vector<double>*> transfer_velocities;
3435 * transfer_velocities.push_back(&transfer_velocity);
3436 * transfer_velocities.push_back(&transfer_velocity_minus_1);
3437 * solution_transfer_velocity.interpolate(transfer_velocities);
3438 * transfer_velocity.update_ghost_values();
3439 * transfer_velocity_minus_1.update_ghost_values();
3440 * solution_transfer_pressure.interpolate(transfer_pressure);
3441 * transfer_pressure.update_ghost_values();
3443 * u_n = transfer_velocity;
3444 * u_n_minus_1 = transfer_velocity_minus_1;
3445 * pres_n = transfer_pressure;
3451 * Interpolate the locally refined solution to a mesh with maximal resolution
3452 * and transfer velocity and pressure.
3459 *
void NavierStokesProjection<dim>::interpolate_max_res(
const unsigned int level) {
3461 * solution_transfer_velocity(dof_handler_velocity);
3462 * std::vector<const LinearAlgebra::distributed::Vector<double>*> velocities;
3463 * velocities.push_back(&u_n);
3464 * velocities.push_back(&u_n_minus_1);
3465 * solution_transfer_velocity.prepare_for_coarsening_and_refinement(velocities);
3468 * solution_transfer_pressure(dof_handler_pressure);
3469 * solution_transfer_pressure.prepare_for_coarsening_and_refinement(pres_n);
3472 *
if(cell->is_locally_owned())
3473 * cell->set_refine_flag();
3480 * transfer_pressure;
3482 * transfer_velocity.
reinit(u_n);
3483 * transfer_velocity.zero_out_ghost_values();
3484 * transfer_velocity_minus_1.reinit(u_n_minus_1);
3485 * transfer_velocity_minus_1.zero_out_ghost_values();
3487 * transfer_pressure.reinit(pres_n);
3488 * transfer_pressure.zero_out_ghost_values();
3490 * std::vector<LinearAlgebra::distributed::Vector<double>*> transfer_velocities;
3492 * transfer_velocities.push_back(&transfer_velocity);
3493 * transfer_velocities.push_back(&transfer_velocity_minus_1);
3494 * solution_transfer_velocity.interpolate(transfer_velocities);
3495 * transfer_velocity.update_ghost_values();
3496 * transfer_velocity_minus_1.update_ghost_values();
3498 * solution_transfer_pressure.interpolate(transfer_pressure);
3499 * transfer_pressure.update_ghost_values();
3501 * u_n = transfer_velocity;
3502 * u_n_minus_1 = transfer_velocity_minus_1;
3503 * pres_n = transfer_pressure;
3509 * Save maximum resolution to a mesh adapted.
3516 *
void NavierStokesProjection<dim>::save_max_res() {
3518 *
GridGenerator::plate_with_a_hole(triangulation_tmp, 0.5, 1.0, 1.0, 1.1, 1.0, 19.0,
Point<2>(2.0, 2.0), 0, 1, 1.0, 2,
true);
3519 * triangulation_tmp.refine_global(
triangulation.n_global_levels() - 1);
3523 * dof_handler_velocity_tmp.distribute_dofs(fe_velocity);
3524 * dof_handler_pressure_tmp.distribute_dofs(fe_pressure);
3528 * u_n_tmp.
reinit(dof_handler_velocity_tmp.n_dofs());
3529 * pres_n_tmp.reinit(dof_handler_pressure_tmp.n_dofs());
3532 * std::vector<std::string> velocity_names(dim,
"v");
3533 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
3536 * u_n_tmp.update_ghost_values();
3537 * data_out.add_data_vector(dof_handler_velocity_tmp, u_n_tmp, velocity_names, component_interpretation_velocity);
3539 * pres_n_tmp.update_ghost_values();
3541 * PostprocessorVorticity<dim> postprocessor;
3542 * data_out.add_data_vector(dof_handler_velocity_tmp, u_n_tmp, postprocessor);
3545 *
const std::string output =
"./" + saving_dir +
"/solution_max_res_end.vtu";
3546 * data_out.write_vtu_in_parallel(output, MPI_COMM_WORLD);
3554 * @sect{ <code>NavierStokesProjection::run</code> }
3558 * This is the time marching function, which starting at <code>t_0</code>
3559 * advances in time
using the projection method with time step <code>dt</code>
3560 * until <code>T</code>.
3564 * Its
second parameter, <code>verbose</code> indicates whether the function
3565 * should output information what it is doing at any given moment:
3573 *
void NavierStokesProjection<dim>::run(
const bool verbose,
const unsigned int output_interval) {
3576 * output_results(1);
3577 *
double time = t_0 + dt;
3578 *
unsigned int n = 1;
3579 *
while(
std::abs(T - time) > 1e-10) {
3582 * pcout <<
"Step = " << n <<
" Time = " << time << std::endl;
3585 * TR_BDF2_stage = 1;
3586 * navier_stokes_matrix.set_TR_BDF2_stage(TR_BDF2_stage);
3588 * mg_matrices[
level].set_TR_BDF2_stage(TR_BDF2_stage);
3590 * verbose_cout <<
" Interpolating the velocity stage 1" << std::endl;
3591 * interpolate_velocity();
3593 * verbose_cout <<
" Diffusion Step stage 1 " << std::endl;
3596 * verbose_cout <<
" Projection Step stage 1" << std::endl;
3598 * u_tmp.equ(gamma*dt, u_tmp);
3600 * projection_step();
3602 * verbose_cout <<
" Updating the Velocity stage 1" << std::endl;
3603 * u_n_gamma.equ(1.0, u_star);
3605 * grad_pres_int.equ(1.0, u_tmp);
3606 * u_tmp.equ(-gamma*dt, u_tmp);
3607 * u_n_gamma += u_tmp;
3608 * u_n_minus_1 = u_n;
3611 * TR_BDF2_stage = 2;
3613 * mg_matrices[
level].set_TR_BDF2_stage(TR_BDF2_stage);
3614 * navier_stokes_matrix.set_TR_BDF2_stage(TR_BDF2_stage);
3616 * verbose_cout <<
" Interpolating the velocity stage 2" << std::endl;
3617 * interpolate_velocity();
3619 * verbose_cout <<
" Diffusion Step stage 2 " << std::endl;
3622 * verbose_cout <<
" Projection Step stage 2" << std::endl;
3623 * u_tmp.equ((1.0 - gamma)*dt, grad_pres_int);
3625 * projection_step();
3627 * verbose_cout <<
" Updating the Velocity stage 2" << std::endl;
3628 * u_n.equ(1.0, u_star);
3630 * u_tmp.equ((gamma - 1.0)*dt, u_tmp);
3633 *
const double max_vel = get_maximal_velocity();
3634 * pcout<<
"Maximal velocity = " << max_vel << std::endl;
3636 * pcout <<
"CFL = " << dt*max_vel*(EquationData::degree_p + 1)*
3638 * compute_lift_and_drag();
3639 *
if(n % output_interval == 0) {
3640 * verbose_cout <<
"Plotting Solution final" << std::endl;
3641 * output_results(n);
3644 *
if(T - time < dt && T - time > 1e-10) {
3646 * navier_stokes_matrix.set_dt(dt);
3648 * mg_matrices[
level].set_dt(dt);
3651 *
if(refinement_iterations > 0 && n % refinement_iterations == 0) {
3652 * verbose_cout <<
"Refining mesh" << std::endl;
3656 *
if(n % output_interval != 0) {
3657 * verbose_cout <<
"Plotting Solution final" << std::endl;
3658 * output_results(n);
3660 *
if(refinement_iterations > 0) {
3661 *
for(
unsigned int lev = 0; lev <
triangulation.n_global_levels() - 1; ++ lev)
3662 * interpolate_max_res(lev);
3674 * @sect{ The main function }
3678 * The main function looks very much like in all the other tutorial programs. We
first initialize MPI,
3679 * we initialize the
class 'NavierStokesProjection' with the dimension as template parameter and then
3680 * let the method
'run' do the job.
3686 *
int main(
int argc,
char *argv[]) {
3688 *
using namespace NS_TRBDF2;
3690 * RunTimeParameters::Data_Storage data;
3691 * data.read_data(
"parameter-file.prm");
3698 * NavierStokesProjection<2> test(data);
3699 * test.run(data.verbose, data.output_interval);
3701 *
if(curr_rank == 0)
3702 * std::cout <<
"----------------------------------------------------"
3704 * <<
"Apparently everything went fine!" << std::endl
3705 * <<
"Don't forget to brush your teeth :-)" << std::endl
3710 *
catch(std::exception &exc) {
3711 * std::cerr << std::endl
3713 * <<
"----------------------------------------------------"
3715 * std::cerr <<
"Exception on processing: " << std::endl
3716 * << exc.what() << std::endl
3717 * <<
"Aborting!" << std::endl
3718 * <<
"----------------------------------------------------"
3723 * std::cerr << std::endl
3725 * <<
"----------------------------------------------------"
3727 * std::cerr <<
"Unknown exception!" << std::endl
3728 * <<
"Aborting!" << std::endl
3729 * <<
"----------------------------------------------------"
3738<a name=
"ann-runtime_parameters.h"></a>
3739<h1>Annotated version of runtime_parameters.h</h1>
3743 * We start by including all the necessary deal.II header files
3749 * #include <deal.II/base/parameter_handler.h>
3755 * @sect{Run time parameters}
3759 * Since our method has several parameters that can be fine-tuned we put them
3760 * into an external file, so that they can be determined at
run-time.
3766 *
namespace RunTimeParameters {
3767 *
using namespace dealii;
3769 *
class Data_Storage {
3773 *
void read_data(
const std::string& filename);
3775 *
double initial_time;
3776 *
double final_time;
3781 *
unsigned int n_refines;
3782 *
unsigned int max_loc_refinements;
3783 *
unsigned int min_loc_refinements;
3787 *
unsigned int max_iterations;
3791 *
unsigned int output_interval;
3795 *
unsigned int refinement_iterations;
3803 * In the constructor of
this class we declare all the parameters in suitable (but arbitrary) subsections.
3809 * Data_Storage::Data_Storage(): initial_time(0.0),
3814 * max_loc_refinements(0),
3815 * min_loc_refinements(0),
3816 * max_iterations(1000),
3819 * output_interval(15),
3820 * refinement_iterations(0) {
3821 * prm.enter_subsection(
"Physical data");
3823 * prm.declare_entry(
"initial_time",
3826 *
" The initial time of the simulation. ");
3827 * prm.declare_entry(
"final_time",
3830 *
" The final time of the simulation. ");
3831 * prm.declare_entry(
"Reynolds",
3834 *
" The Reynolds number. ");
3836 * prm.leave_subsection();
3838 * prm.enter_subsection(
"Time step data");
3840 * prm.declare_entry(
"dt",
3843 *
" The time step size. ");
3845 * prm.leave_subsection();
3847 * prm.enter_subsection(
"Space discretization");
3849 * prm.declare_entry(
"n_of_refines",
3852 *
" The number of cells we want on each direction of the mesh. ");
3853 * prm.declare_entry(
"max_loc_refinements",
3856 *
" The number of maximum local refinements. ");
3857 * prm.declare_entry(
"min_loc_refinements",
3860 *
" The number of minimum local refinements. ");
3862 * prm.leave_subsection();
3864 * prm.enter_subsection(
"Data solve");
3866 * prm.declare_entry(
"max_iterations",
3869 *
" The maximal number of iterations linear solvers must make. ");
3870 * prm.declare_entry(
"eps",
3873 *
" The stopping criterion. ");
3875 * prm.leave_subsection();
3877 * prm.declare_entry(
"refinement_iterations",
3880 *
" This number indicates how often we need to "
3881 *
"refine the mesh");
3883 * prm.declare_entry(
"saving directory",
"SimTest");
3885 * prm.declare_entry(
"verbose",
3888 *
" This indicates whether the output of the solution "
3889 *
"process should be verbose. ");
3891 * prm.declare_entry(
"output_interval",
3894 *
" This indicates between how many time steps we print "
3895 *
"the solution. ");
3900 * We need now a routine to read all declared parameters in the constructor
3906 *
void Data_Storage::read_data(
const std::string& filename) {
3907 * std::ifstream file(filename);
3910 * prm.parse_input(file);
3912 * prm.enter_subsection(
"Physical data");
3914 * initial_time = prm.get_double(
"initial_time");
3915 * final_time = prm.get_double(
"final_time");
3916 * Reynolds = prm.get_double(
"Reynolds");
3918 * prm.leave_subsection();
3920 * prm.enter_subsection(
"Time step data");
3922 * dt = prm.get_double(
"dt");
3924 * prm.leave_subsection();
3926 * prm.enter_subsection(
"Space discretization");
3928 * n_refines = prm.get_integer(
"n_of_refines");
3929 * max_loc_refinements = prm.get_integer(
"max_loc_refinements");
3930 * min_loc_refinements = prm.get_integer(
"min_loc_refinements");
3932 * prm.leave_subsection();
3934 * prm.enter_subsection(
"Data solve");
3936 * max_iterations = prm.get_integer(
"max_iterations");
3937 *
eps = prm.get_double(
"eps");
3939 * prm.leave_subsection();
3941 * dir = prm.get(
"saving directory");
3943 * refinement_iterations = prm.get_integer(
"refinement_iterations");
3945 * verbose = prm.get_bool(
"verbose");
3947 * output_interval = prm.get_integer(
"output_interval");
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
unsigned int depth_console(const unsigned int n)
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&...args)
void build(const DoFHandler< dim, dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
void loop(const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &face_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &boundary_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
types::boundary_id get_boundary_id(const unsigned int face_batch_index) const
void initialize_dof_vector(VectorType &vec, const unsigned int dof_handler_index=0) const
void cell_loop(const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false) const
__global__ void set(Number *val, const Number s, const size_type N)
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
void loop(ITERATOR begin, std_cxx20::type_identity_t< ITERATOR > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
void mesh_loop(const CellIteratorType &begin, const CellIteratorType &end, const CellWorkerFunctionType &cell_worker, const CopierType &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const BoundaryWorkerFunctionType &boundary_worker=BoundaryWorkerFunctionType(), const FaceWorkerFunctionType &face_worker=FaceWorkerFunctionType(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
@ component_is_part_of_vector
void create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
void plate_with_a_hole(Triangulation< dim > &tria, const double inner_radius=0.4, const double outer_radius=1., const double pad_bottom=2., const double pad_top=2., const double pad_left=1., const double pad_right=1., const Point< dim > ¢er=Point< dim >(), const types::manifold_id polar_manifold_id=0, const types::manifold_id tfi_manifold_id=1, const double L=1., const unsigned int n_slices=2, const bool colorize=false)
Rectangular plate with an (offset) cylindrical hole.
@ matrix
Contents is actually a matrix.
@ diagonal
Matrix is diagonal.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
@ construct_multigrid_hierarchy
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
long double gamma(const unsigned int n)
int(&) functions(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
static const unsigned int invalid_unsigned_int
void refine_and_coarsen_fixed_number(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const types::global_cell_index max_n_cells=std::numeric_limits< types::global_cell_index >::max())
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector<::Vector< double > > solution_values
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients
TasksParallelScheme tasks_parallel_scheme
UpdateFlags mapping_update_flags
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
std_cxx20::type_identity< T > identity
const ::Triangulation< dim, spacedim > & tria