8#ifndef OPENFPM_NUMERICS_SRC_FINITEDIFFERENCE_EQ_UNIT_TEST_3D_HPP_
9#define OPENFPM_NUMERICS_SRC_FINITEDIFFERENCE_EQ_UNIT_TEST_3D_HPP_
11#define BOOST_TEST_DYN_LINK
12#include <boost/test/unit_test.hpp>
15#include "Laplacian.hpp"
16#include "FiniteDifference/eq.hpp"
17#include "FiniteDifference/sum.hpp"
18#include "FiniteDifference/mul.hpp"
19#include "Grid/grid_dist_id.hpp"
20#include "Decomposition/CartDecomposition.hpp"
21#include "Vector/Vector.hpp"
22#include "Solvers/umfpack_solver.hpp"
23#include "data_type/aggregate.hpp"
24#include "Solvers/petsc_solver.hpp"
25#include "FiniteDifference/FDScheme.hpp"
27BOOST_AUTO_TEST_SUITE( eq_test_suite_3d )
33 static const unsigned int dims = 3;
35 static const unsigned int nvar = 4;
38 static const bool boundary[];
59 static const unsigned int dims = 3;
61 static const unsigned int nvar = 4;
94 static float val() {
return 1.0;}
97template<
typename solver_type,
typename l
id_nn_3d>
void lid_driven_cavity_3d()
99 #include "Equations/stoke_flow_eq_3d.hpp"
112 long int sz[] = {36,12,12};
114 szu[0] = (size_t)sz[0];
115 szu[1] = (size_t)sz[1];
116 szu[2] = (size_t)sz[2];
121 constexpr int velocity = 0;
122 constexpr int pressure = 1;
135 fd.impose(
ic_eq(),0.0, EQ_4, {0,0,0},{sz[0]-2,sz[1]-2,sz[2]-2},
true);
136 fd.impose(
Prs(), 0.0, EQ_4, {0,0,0},{0,0,0});
137 fd.impose(
vx_eq(),0.0, EQ_1, {1,0},{sz[0]-2,sz[1]-2,sz[2]-2});
138 fd.impose(
vy_eq(),0.0, EQ_2, {0,1},{sz[0]-2,sz[1]-2,sz[2]-2});
139 fd.impose(
vz_eq(),0.0, EQ_3, {0,0,1},{sz[0]-2,sz[1]-2,sz[2]-2});
143 fd.impose(
v_x(),0.0, EQ_1, {0,0,0}, {0,sz[1]-2,sz[2]-2});
144 fd.impose(
v_x(),0.0, EQ_1, {sz[0]-1,0,0},{sz[0]-1,sz[1]-2,sz[2]-2});
147 fd.impose(
avg_y_vx_f(),0.0, EQ_1, {0,-1,0}, {sz[0]-1,-1,sz[2]-2});
148 fd.impose(
avg_y_vx(),0.0, EQ_1, {0,sz[1]-1,0},{sz[0]-1,sz[1]-1,sz[2]-2});
151 fd.impose(
avg_z_vx_f(),0.0, EQ_1, {0,-1,-1}, {sz[0]-1,sz[1]-1,-1});
152 fd.impose(
avg_z_vx(),0.0, EQ_1, {0,-1,sz[2]-1},{sz[0]-1,sz[1]-1,sz[2]-1});
156 fd.impose(
avg_x_vy_f(),0.0, EQ_2, {-1,0,0}, {-1,sz[1]-1,sz[2]-2});
157 fd.impose(
avg_x_vy(),1.0, EQ_2, {sz[0]-1,0,0},{sz[0]-1,sz[1]-1,sz[2]-2});
160 fd.impose(
v_y(), 0.0, EQ_2, {0,0,0}, {sz[0]-2,0,sz[2]-2});
161 fd.impose(
v_y(), 0.0, EQ_2, {0,sz[1]-1,0},{sz[0]-2,sz[1]-1,sz[2]-2});
164 fd.impose(
avg_z_vy(),0.0, EQ_2, {-1,0,sz[2]-1}, {sz[0]-1,sz[1]-1,sz[2]-1});
165 fd.impose(
avg_z_vy_f(),0.0, EQ_2, {-1,0,-1}, {sz[0]-1,sz[1]-1,-1});
169 fd.impose(
avg_x_vz_f(),0.0, EQ_3, {-1,0,0}, {-1,sz[1]-2,sz[2]-1});
170 fd.impose(
avg_x_vz(),1.0, EQ_3, {sz[0]-1,0,0},{sz[0]-1,sz[1]-2,sz[2]-1});
173 fd.impose(
avg_y_vz(),0.0, EQ_3, {-1,sz[1]-1,0},{sz[0]-1,sz[1]-1,sz[2]-1});
174 fd.impose(
avg_y_vz_f(),0.0, EQ_3, {-1,-1,0}, {sz[0]-1,-1,sz[2]-1});
177 fd.impose(
v_z(),0.0, EQ_3, {0,0,0}, {sz[0]-2,sz[1]-2,0});
178 fd.impose(
v_z(),0.0, EQ_3, {0,0,sz[2]-1},{sz[0]-2,sz[1]-2,sz[2]-1});
183 fd.impose(
Prs(), 0.0, EQ_4, {-1,-1,-1},{-1,sz[1]-1,sz[2]-1});
184 fd.impose(
Prs(), 0.0, EQ_4, {sz[0]-1,-1,-1},{sz[0]-1,sz[1]-1,sz[2]-1});
187 fd.impose(
Prs(), 0.0, EQ_4, {0,sz[1]-1,-1}, {sz[0]-2,sz[1]-1,sz[2]-1});
188 fd.impose(
Prs(), 0.0, EQ_4, {0,-1 ,-1}, {sz[0]-2,-1, sz[2]-1});
191 fd.impose(
Prs(), 0.0, EQ_4, {0,0,sz[2]-1}, {sz[0]-2,sz[1]-2,sz[2]-1});
192 fd.impose(
Prs(), 0.0, EQ_4, {0,0,-1}, {sz[0]-2,sz[1]-2,-1});
195 fd.impose(
v_x(), 0.0, EQ_1, {-1,-1,-1},{-1,sz[1]-1,sz[2]-1});
196 fd.impose(
v_y(), 0.0, EQ_2, {-1,-1,-1},{sz[0]-1,-1,sz[2]-1});
197 fd.impose(
v_z(), 0.0, EQ_3, {-1,-1,-1},{sz[0]-1,sz[1]-1,-1});
200 auto x_ = solver.try_solve(fd.getA(),fd.getB());
203 fd.template copy<velocity,pressure>(x_,{0,0},{sz[0]-1,sz[1]-1,sz[2]-1},g_dist);
205 std::string s = std::string(demangle(
typeid(solver_type).name()));
208 g_dist.write(s +
"lid_driven_cavity_3d_p" + std::to_string(v_cl.
getProcessingUnits()) +
"_grid");
210#if !(defined(SE_CLASS3) || defined(TEST_COVERAGE_MODE))
214 g_dist2.load(
"test/lid_driven_cavity_3d_reference.hdf5");
216 auto it2 = g_dist2.getDomainIterator();
223 test &= fabs(g_dist2.template getProp<velocity>(p)[0] - g_dist.template getProp<velocity>(p)[0]) < 3.5e-5;
224 test &= fabs(g_dist2.template getProp<velocity>(p)[1] - g_dist.template getProp<velocity>(p)[1]) < 3.5e-5;
225 test &= fabs(g_dist2.template getProp<velocity>(p)[2] - g_dist.template getProp<velocity>(p)[2]) < 3.5e-5;
227 test &= fabs(g_dist2.template getProp<pressure>(p) - g_dist.template getProp<pressure>(p)) < 3.0e-4;
231 std::cout << g_dist2.template getProp<velocity>(p)[0] <<
" " << g_dist.template getProp<velocity>(p)[0] << std::endl;
232 std::cout << g_dist2.template getProp<velocity>(p)[1] <<
" " << g_dist.template getProp<velocity>(p)[1] << std::endl;
233 std::cout << g_dist2.template getProp<velocity>(p)[2] <<
" " << g_dist.template getProp<velocity>(p)[2] << std::endl;
235 std::cout << g_dist2.template getProp<pressure>(p) <<
" " << g_dist.template getProp<pressure>(p) << std::endl;
243 BOOST_REQUIRE_EQUAL(test,
true);
251BOOST_AUTO_TEST_CASE(lid_driven_cavity)
253#if defined(HAVE_EIGEN) && defined(HAVE_SUITESPARSE)
261BOOST_AUTO_TEST_SUITE_END()
This class represent an N-dimensional box.
This class decompose a space into sub-sub-domains and distribute them across processors.
Class that contain Padding information on each direction positive and Negative direction.
Sparse Matrix implementation.
size_t getProcessingUnits()
Get the total number of processors.
Implementation of VCluster class.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
This is a distributed grid.
[Definition of the system]
static float val()
therutn the value of the constant
void const_field
define that eta is a constant field
Specify the general caratteristic of system to solve.
SparseMatrix< double, int, EIGEN_BASE > SparseMatrix_type
type of SparseMatrix for the linear solver
static const bool boundary[]
boundary at X and Y
float stype
type of space float, double, ...
Vector< double > Vector_type
type of Vector for the linear solver
grid_dist_id< 3, float, aggregate< float[3], float >, CartDecomposition< 3, float > > b_grid
type of base grid
static const unsigned int nvar
number of fields in the system
grid_dist_id< 3, float, aggregate< float[3], float >, CartDecomposition< 3, float > > b_grid
type of base grid
static const bool boundary[]
boundary at X and Y
float stype
type of space float, double, ...
static const unsigned int dims
dimensionaly of the equation ( 3D problem ...)
Vector< double, PETSC_BASE > Vector_type
type of Vector for the linear solver
SparseMatrix< double, int, PETSC_BASE > SparseMatrix_type
type of SparseMatrix for the linear solver
It model an expression expr1 + ... exprn.