OpenLB 1.7
Loading...
Searching...
No Matches
wallFunctionBoundaryPostProcessors3D.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2018 Marc Haussmann
4 * E-mail contact: info@openlb.net
5 * The most recent release of OpenLB can be downloaded at
6 * <http://www.openlb.net/>
7 *
8 * This program is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * as published by the Free Software Foundation; either version 2
11 * of the License, or (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public
19 * License along with this program; if not, write to the Free
20 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
21 * Boston, MA 02110-1301, USA.
22*/
23
24#ifndef WALLFUNCTION_BOUNDARY_POST_PROCESSORS_3D_HH
25#define WALLFUNCTION_BOUNDARY_POST_PROCESSORS_3D_HH
26
29#include "dynamics/lbm.h"
30#include "core/util.h"
33
34namespace olb {
35
36template <typename T, typename S>
37Musker<T,S>::Musker(T nu, T y, T rho) : AnalyticalF<1,T,S>(1), _nu(nu), _y(y),_rho(rho)
38{
39 this->getName() = "Musker";
40}
41
42template <typename T, typename S>
43bool Musker<T,S>::operator()(T output[], const S tau_w[])
44{
45 T y_plus = _y*util::sqrt(tau_w[0]/_rho)/_nu;
46
47 T a = 5.424;
48 T b = 0.119760479041916168;
49 T c = 0.488023952095808383;
50 T d = 0.434;
51 T e = 3.50727901936264842;
52
53 output[0] = util::sqrt(tau_w[0]/_rho)*(a*util::atan(b*y_plus - c) +
54 d*util::log(util::pow(y_plus+10.6, 9.6)/util::pow(util::pow(y_plus, 2) - 8.15*y_plus + 86, 2)) - e);
55
56 // Condition for the sub-viscous layer : TODO MARC H
57 if (output[0] < 0) {
58 output[0] = y_plus * util::sqrt(tau_w[0]/_rho);
59 }
60
61 return true;
62}
63
64template <typename T, typename S>
65PowerLawProfile<T,S>::PowerLawProfile(T nu, T u2, T y2, T y1, T rho) : AnalyticalF<1,T,S>(2), _nu(nu), _u2(u2), _y2(y2), _y1(y1), _rho(rho)
66{
67 this->getName() = "PowerLawProfile";
68}
69
70template <typename T, typename S>
71bool PowerLawProfile<T,S>::operator()(T output[], const S input[])
72{
73 T tau_w = 0.0246384 * _rho * util::pow(_nu, 0.25) * util::pow(_u2, 1.75) / util::pow(_y2, 0.25);
74 T u_tau = util::sqrt(tau_w/_rho);
75 T y_plus = _y1 * u_tau / _nu;
76
77 if (y_plus > 30.0) {
78 output[0] = tau_w;
79 output[1] = u_tau * 8.3 * util::pow(y_plus, 1./7.);
80 }
81 else if (y_plus < 30.0 && y_plus > 5.) {
82 output[0] = tau_w;
83 output[1] = u_tau * (-2.81742 + 4.85723 * util::log(y_plus));
84 }
85 else {
86 output[0] = 2.*_u2 * _rho * _nu / _y2;
87 if (output[0] < 0.) {
88 output[0]*=-1.;
89 }
90 u_tau = util::sqrt(output[0]/_rho);
91 y_plus = _y1 * u_tau / _nu;
92 output[1] = u_tau * y_plus;
93 }
94 return true;
95}
96
97
99template<typename T, typename DESCRIPTOR>
101 BlockGeometry<T,3>& blockGeometryStructure,
102 std::vector<int> discreteNormal, std::vector<int> missingIndices,
104 IndicatorF3D<T>* geoIndicator)
105 : x0(x0_), x1(x1_), y0(y0_), y1(y1_), z0(z0_), z1(z1_),
106 _blockGeometryStructure(blockGeometryStructure),
107 _discreteNormal(discreteNormal), _missingIndices(missingIndices),
108 _converter(converter), _wallFunctionParam(wallFunctionParam)
109{
110 this->getName() = "WallFunctionBoundaryProcessor3D";
111 // Define Direction and orientatation
112 discreteNormalX = _discreteNormal[0];
113 discreteNormalY = _discreteNormal[1];
114 discreteNormalZ = _discreteNormal[2];
115
116 int Type_BC = discreteNormalX*discreteNormalX + discreteNormalY*discreteNormalY + discreteNormalZ*discreteNormalZ;
117 T normal_norm = util::sqrt(Type_BC); // l2 norm : magnitude of the vector
118 if (Type_BC == 1) { // Straight plane
119 if (discreteNormalX != 0) {
120 orientation = discreteNormalX;
121 direction = 0;
122 }
123 else if (discreteNormalY != 0) {
124 orientation = discreteNormalY;
125 direction = 1;
126 }
127 else if (discreteNormalZ != 0) {
128 orientation = discreteNormalZ;
129 direction = 2;
130 }
131 }
132 else if (Type_BC == 2) { // Edge
133 orientation = 0;
134 direction = 0;
135 }
136 else if (Type_BC == 3) { // Corner
137 orientation = 0;
138 direction = 0;
139 }
140
141 if (geoIndicator == NULL) {
143 y_1 = 0.5; // Default half-way Bounce-Back distance - [m] DESCRIPTOR units
144 if (_wallFunctionParam.latticeWalldistance > 0.) {
145 y_1 = _wallFunctionParam.latticeWalldistance; // [m] DESCRIPTOR units
146 }
147 y_1 *= normal_norm; // DESCRIPTOR units
148
149 // Define the unit normal vector
150 unit_normal[0] = discreteNormalX / normal_norm;
151 unit_normal[1] = discreteNormalY / normal_norm;
152 unit_normal[2] = discreteNormalZ / normal_norm;
153 }
154 else {
155 calculateWallDistances(geoIndicator);
156 if (y_1 > 2. || std::isnan(unit_normal[0]) || std::isnan(unit_normal[1]) || std::isnan(unit_normal[2])) {
157 y_1 = 0.5; // Default half-way Bounce-Back distance - [m] DESCRIPTOR units
158 if (_wallFunctionParam.latticeWalldistance > 0.) {
159 y_1 = _wallFunctionParam.latticeWalldistance; // [m] DESCRIPTOR units
160 }
161 y_1 *= normal_norm; // DESCRIPTOR units
162
163 // Define the unit normal vector
164 unit_normal[0] = discreteNormalX / normal_norm;
165 unit_normal[1] = discreteNormalY / normal_norm;
166 unit_normal[2] = discreteNormalZ / normal_norm;
167 }
168 }
169 y_2 = y_1 + normal_norm; // DESCRIPTOR units
170
171 // Computation of Rho - Zou an He - Speed up
172 getIndices(direction, 0, onWallIndices);
173 getIndices(direction, orientation, normalIndices);
174
175 // Vector needed for Pi Computation - FNeq Bounce-Back - Speep up
176 // Malaspinas condition c*n < zero
177 // Definition of directions pointing towards the fluid
178 getIndices(direction, -orientation, normalInwardsIndices);
179}
180
181template<typename T, typename DESCRIPTOR>
182void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::processSubDomain(BlockLattice<T,DESCRIPTOR>& blockLattice, int x0_, int x1_, int y0_, int y1_, int z0_, int z1_)
183{
184 int newX0, newX1, newY0, newY1, newZ0, newZ1;
185 if ( util::intersect (
186 x0, x1, y0, y1, z0, z1,
187 x0_, x1_, y0_, y1_, z0_, z1_,
188 newX0, newX1, newY0, newY1, newZ0, newZ1 ) ) {
189
190 int iX;
191#ifdef PARALLEL_MODE_OMP
192 #pragma omp parallel for
193#endif
194 for (iX=newX0; iX<=newX1; ++iX) {
195 for (int iY=newY0; iY<=newY1; ++iY) {
196 for (int iZ=newZ0; iZ<=newZ1; ++iZ) {
197 ComputeWallFunction(blockLattice, iX,iY,iZ);
198 }
199 }
200 }
201 }
202}
203
204template<typename T, typename DESCRIPTOR>
205void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::getIndices(int index, int value, std::vector<int>& indices)
206{
207 for (int iVel=0; iVel<DESCRIPTOR::q; ++iVel) {
208 if (descriptors::c<DESCRIPTOR>(iVel,index)==value) {
209 indices.push_back(iVel);
210 }
211 }
212}
213
214template<typename T, typename DESCRIPTOR>
215void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::calculateWallDistances(IndicatorF3D<T>* indicator)
216{
217 using S = BaseType<T>;
218
219 T physR[3];
220 int iX = x0;
221 int iY = y1;
222 int iZ = z1;
223 T scaling = _converter.getConversionFactorLength() * 0.1;
224 _blockGeometryStructure.getPhysR(physR,{iX, iY, iZ});
225 Vector<T,3> origin(physR[0],physR[1],physR[2]);
226 Vector<T,3> normal(0.,0.,0.);
227 T distance = 0.;
228 Vector<T,3> direction(0.,0.,0.);
229 int smallestDistance_i = 0;
230 T smallestDistance = 0.;
231 bool firstHit = true;
232 origin[0] = physR[0];
233 origin[1] = physR[1];
234 origin[2] = physR[2];
235 int discreteDirection[6][3];
236 for (int i=0; i < 6; i++) {
237 for (int j=0; j < 3; j++) {
238 discreteDirection[i][j] = 0;
239 }
240 }
241 discreteDirection[0][0] = 1;
242 discreteDirection[1][0] = -1;
243 discreteDirection[2][1] = 1;
244 discreteDirection[3][1] = -1;
245 discreteDirection[4][2] = 1;
246 discreteDirection[5][2] = -1;
247 for (int i=0; i < 6; i++) {
248 direction[0] = discreteDirection[i][0] * scaling;
249 direction[1] = discreteDirection[i][1] * scaling;
250 direction[2] = discreteDirection[i][2] * scaling;
251 if (indicator->distance(distance, origin, direction)) {
252 if (firstHit) {
253 smallestDistance = distance;
254 smallestDistance_i = i;
255 firstHit = false;
256 }
257 else {
258 if (distance < smallestDistance) {
259 smallestDistance = distance;
260 smallestDistance_i = i;
261 }
262 }
263 }
264 }
265
266 direction[0] = discreteDirection[smallestDistance_i][0] * scaling;
267 direction[1] = discreteDirection[smallestDistance_i][1] * scaling;
268 direction[2] = discreteDirection[smallestDistance_i][2] * scaling;
269
270 Vector<T,3> direction2(direction);
271 Vector<T,3> direction3(direction);
272 if (smallestDistance_i == 0 || smallestDistance_i == 1 ) {
273 direction2[1] = direction2[0] * scaling;
274 direction3[2] = direction3[0] * scaling;
275 }
276 else if (smallestDistance_i == 2 || smallestDistance_i == 3 ) {
277 direction2[0] = direction2[1] * scaling;
278 direction3[2] = direction3[1] * scaling;
279 }
280 else {
281 direction2[0] = direction2[2] * scaling;
282 direction3[1] = direction3[2] * scaling;
283 }
284
285 Vector<S,3> directionN = direction*(1/norm(direction));
286 Vector<S,3> POS(origin + smallestDistance*directionN); //Point on Surface
287
288 indicator->distance(distance, origin, direction2);
289 Vector<S,3> direction2N = direction2*(1/norm(direction2));
290 Vector<S,3> POS2(origin + distance*direction2N); //Point on Surface
291
292 indicator->distance(distance, origin, direction3);
293
294 Vector<S,3> direction3N = direction3*(1/norm(direction3));
295 Vector<S,3> POS3(origin + distance*direction3N); //Point on Surface
296
297 Vector<S,3> vec1 (POS - POS2);
298 Vector<S,3> vec2 (POS - POS3);
299
300 normal[0] = -(vec1[1]*vec2[2] - vec1[2]*vec2[1]);
301 normal[1] = -(vec1[2]*vec2[0] - vec1[0]*vec2[2]);
302 normal[2] = -(vec1[0]*vec2[1] - vec1[1]*vec2[0]);
303
304 T normalMagnitude = util::sqrt(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
305 normal[0] /= normalMagnitude;
306 normal[1] /= normalMagnitude;
307 normal[2] /= normalMagnitude;
308
309 unit_normal[0] = normal[0];
310 unit_normal[1] = normal[0];
311 unit_normal[2] = normal[0];
312
313 direction[0] = normal[0] * scaling;
314 direction[1] = normal[1] * scaling;
315 direction[2] = normal[2] * scaling;
316
317 indicator->distance(distance, origin, direction);
318 y_1 = distance / _converter.getConversionFactorLength();
319}
320
321template<typename T, typename DESCRIPTOR>
322void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::VelGradFromSecondOrderFD(bool NormalGradient, T Vel_BC[DESCRIPTOR::d], T Vel_1[DESCRIPTOR::d], T Vel_2[DESCRIPTOR::d], T VelGrad[DESCRIPTOR::d])
323{
324 if (NormalGradient) {
325 if (orientation == 1) {
326 for (int Dim=0; Dim<DESCRIPTOR::d; Dim++) {
327 VelGrad[Dim] = fd::BSGradient(Vel_BC[Dim], Vel_1[Dim], Vel_2[Dim]); // Backward second-order velocity gradient
328 }
329 }
330 else { // orientation == -1
331 for (int Dim=0; Dim<DESCRIPTOR::d; Dim++) {
332 VelGrad[Dim] = fd::FSGradient(Vel_BC[Dim], Vel_1[Dim], Vel_2[Dim]); // Forward second-order velocity gradient
333 }
334 }
335 }
336 else {
337 if (orientation == 1) {
338 for (int Dim=0; Dim<DESCRIPTOR::d; Dim++) {
339 VelGrad[Dim] = fd::centralGradient(Vel_1[Dim], Vel_2[Dim]); // central second-order velocity gradient: centralGradient( y(x+1), y(x-1) )
340 }
341 }
342 else { // orientation == -1
343 for (int Dim=0; Dim<DESCRIPTOR::d; Dim++) {
344 VelGrad[Dim] = fd::centralGradient(Vel_2[Dim], Vel_1[Dim]); // central second-order velocity gradient: centralGradient( y(x+1), y(x-1) )
345 }
346 }
347 }
348}
349
350template<typename T, typename DESCRIPTOR>
351void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeVelocityGradientTensor(T u_bc[DESCRIPTOR::d], T u_x1[DESCRIPTOR::d], T u_x2[DESCRIPTOR::d], T u_y1[DESCRIPTOR::d], T u_y2[DESCRIPTOR::d], T u_z1[DESCRIPTOR::d], T u_z2[DESCRIPTOR::d], T VelGrad[DESCRIPTOR::d][DESCRIPTOR::d])
352{
353 T dx_u[DESCRIPTOR::d], dy_u[DESCRIPTOR::d], dz_u[DESCRIPTOR::d];
354 VelGradFromSecondOrderFD(direction == 0, u_bc, u_x1, u_x2, dx_u);
355 VelGradFromSecondOrderFD(direction == 1, u_bc, u_y1, u_y2, dy_u);
356 VelGradFromSecondOrderFD(direction == 2, u_bc, u_z1, u_z2, dz_u);
357
358 // du/dx du/dy du/dz
359 // dv/dx dv/dy dv/dz
360 // dw/dx dw/dy dw/dz
361 for (int Dim = 0; Dim < DESCRIPTOR::d; Dim++) {
362 VelGrad[Dim][0] = dx_u[Dim];
363 VelGrad[Dim][1] = dy_u[Dim];
364 VelGrad[Dim][2] = dz_u[Dim];
365 }
366}
367
368template<typename T, typename DESCRIPTOR>
369void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeNeighborsU(BlockLattice<T,DESCRIPTOR>& blockLattice, int x, int y, int z,
370 T u_x1[DESCRIPTOR::d], T u_x2[DESCRIPTOR::d], T u_y1[DESCRIPTOR::d], T u_y2[DESCRIPTOR::d], T u_z1[DESCRIPTOR::d], T u_z2[DESCRIPTOR::d])
371{
372 using namespace olb::util::tensorIndices3D;
373 // Computation of local external velocity around lattice node (x,y,z)
374 if (direction == 0) {
375 blockLattice.get(x - discreteNormalX, y, z).computeU(u_x1);
376 blockLattice.get(x - 2*discreteNormalX, y, z).computeU(u_x2);
377 }
378 else {
379 ComputeUWall(blockLattice, x + discreteNormalY + discreteNormalZ, y, z, u_x1);
380 ComputeUWall(blockLattice, x - discreteNormalY - discreteNormalZ, y, z, u_x2);
381 }
382 if (direction == 1) {
383 blockLattice.get(x, y - discreteNormalY, z).computeU(u_y1);
384 blockLattice.get(x, y - 2*discreteNormalY, z).computeU(u_y2);
385 }
386 else {
387 ComputeUWall(blockLattice, x, y + discreteNormalX + discreteNormalZ, z, u_y1);
388 ComputeUWall(blockLattice, x, y - discreteNormalX - discreteNormalZ, z, u_y2);
389 }
390 if (direction == 2) {
391 blockLattice.get(x, y, z - discreteNormalZ).computeU(u_z1);
392 blockLattice.get(x, y, z - 2*discreteNormalZ).computeU(u_z2);
393 }
394 else {
395 ComputeUWall(blockLattice, x, y, z + discreteNormalX + discreteNormalY, u_z1);
396 ComputeUWall(blockLattice, x, y, z - discreteNormalX - discreteNormalY, u_z2);
397 }
398}
399
400template<typename T, typename DESCRIPTOR>
401void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeVelocityGradient(BlockLattice<T,DESCRIPTOR>& blockLattice, int x, int y, int z, T u_bc[DESCRIPTOR::d], T VelGrad[DESCRIPTOR::d][DESCRIPTOR::d])
402{
403 // Computation of neighbor velocity around lattice node (x,y,z)
404 T u_x1[3], u_x2[3], u_y1[3], u_y2[3], u_z1[3], u_z2[3];
405 computeNeighborsU(blockLattice, x, y, z, u_x1, u_x2, u_y1, u_y2, u_z1, u_z2);
406 // Computation of Velocity Gradient Tensor
407 computeVelocityGradientTensor(u_bc, u_x1, u_x2, u_y1, u_y2, u_z1, u_z2, VelGrad);
408}
409
410template<typename T, typename DESCRIPTOR>
412 T u_x1[DESCRIPTOR::d], T u_x2[DESCRIPTOR::d], T u_y1[DESCRIPTOR::d], T u_y2[DESCRIPTOR::d], T u_z1[DESCRIPTOR::d], T u_z2[DESCRIPTOR::d],
413 T& rho_x1, T& rho_x2, T& rho_y1, T& rho_y2, T& rho_z1, T& rho_z2)
414{
415 using namespace olb::util::tensorIndices3D;
416 // Computation of local external velocity around lattice node (x,y,z)
417 if (direction == 0) {
418 rho_x1 = blockLattice.get(x - discreteNormalX, y, z).computeRho();
419 rho_x2 = blockLattice.get(x - 2*discreteNormalX, y, z).computeRho();
420 }
421 else {
422 ComputeRhoWall(blockLattice, blockLattice.get(x + discreteNormalY + discreteNormalZ, y, z), x + discreteNormalY + discreteNormalZ, y, z, u_x1, rho_x1);
423 ComputeRhoWall(blockLattice, blockLattice.get(x - discreteNormalY - discreteNormalZ, y, z), x - discreteNormalY - discreteNormalZ, y, z, u_x2, rho_x2);
424 }
425 if (direction == 1) {
426 rho_y1 = blockLattice.get(x, y - discreteNormalY, z).computeRho();
427 rho_y2 = blockLattice.get(x, y - 2*discreteNormalY, z).computeRho();
428 }
429 else {
430 ComputeRhoWall(blockLattice, blockLattice.get(x, y - discreteNormalX - discreteNormalZ, z), x, y - discreteNormalX - discreteNormalZ, z, u_y1, rho_y1);
431 ComputeRhoWall(blockLattice, blockLattice.get(x, y - discreteNormalX - discreteNormalZ, z), x, y - discreteNormalX - discreteNormalZ, z, u_y2, rho_y2);
432 }
433 if (direction == 2) {
434 rho_z1 = blockLattice.get(x, y, z - discreteNormalZ).computeRho();
435 rho_z2 = blockLattice.get(x, y, z - 2*discreteNormalZ).computeRho();
436 }
437 else {
438 ComputeRhoWall(blockLattice, blockLattice.get(x, y, z + discreteNormalX + discreteNormalY), x, y, z + discreteNormalX + discreteNormalY, u_z1, rho_z1);
439 ComputeRhoWall(blockLattice, blockLattice.get(x, y, z - discreteNormalX - discreteNormalY), x, y, z - discreteNormalX - discreteNormalY, u_z2, rho_z2);
440 }
441}
442
443template<typename T, typename DESCRIPTOR>
445 T u_x1[DESCRIPTOR::d], T u_x2[DESCRIPTOR::d], T u_y1[DESCRIPTOR::d], T u_y2[DESCRIPTOR::d], T u_z1[DESCRIPTOR::d], T u_z2[DESCRIPTOR::d],
446 T& rho_x1, T& rho_x2, T& rho_y1, T& rho_y2, T& rho_z1, T& rho_z2)
447{
448 // Computation of neighbor velocity around lattice node (x,y,z)
449 computeNeighborsU(blockLattice, x, y, z, u_x1, u_x2, u_y1, u_y2, u_z1, u_z2);
450 // Computation of neighbor density around lattice node (x,y,z)
451 computeNeighborsRho(blockLattice, x, y, z, u_x1, u_x2, u_y1, u_y2, u_z1, u_z2, rho_x1, rho_x2, rho_y1, rho_y2, rho_z1, rho_z2);
452}
453
454
456template<typename T, typename DESCRIPTOR>
457void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::ComputeUWall(BlockLattice<T,DESCRIPTOR>& blockLattice, int x, int y, int z, T u[DESCRIPTOR::d])
458{
460 using namespace olb::util::tensorIndices3D;
461 Cell<T,DESCRIPTOR> cell = blockLattice.get(x,y,z);
462
464 T u_2[DESCRIPTOR::d]; // Velocity to the neighbord lattice in the normal inward direction
465 T rho_2; // Density in the neighbor lattice in the inwards direction
466 blockLattice.get(x - discreteNormalX, y - discreteNormalY, z - discreteNormalZ).computeRhoU(rho_2,u_2);
467
469 T u_2_dot_unit_normal = u_2[0] * unit_normal[0] +
470 u_2[1] * unit_normal[1] +
471 u_2[2] * unit_normal[2]; //scalar Product of the u2 and uniy normal
472 T u_2_parallel[3];
473 u_2_parallel[0] = u_2[0] - (u_2_dot_unit_normal * unit_normal[0]);
474 u_2_parallel[1] = u_2[1] - (u_2_dot_unit_normal * unit_normal[1]);
475 u_2_parallel[2] = u_2[2] - (u_2_dot_unit_normal * unit_normal[2]);
476
477 T u_2_parallel_norm = util::sqrt(u_2_parallel[0] * u_2_parallel[0] +
478 u_2_parallel[1] * u_2_parallel[1] +
479 u_2_parallel[2] * u_2_parallel[2]); // l2 norm : magnitude of the vector
480
481 T e_x_loc[3] = {0., 0., 0.}; // Streamwise direction
482 if (u_2_parallel_norm != 0) {
483 T inv_u_2_parallel_norm = (1. /u_2_parallel_norm);
484 e_x_loc[0] = u_2_parallel[0] * inv_u_2_parallel_norm;
485 e_x_loc[1] = u_2_parallel[1] * inv_u_2_parallel_norm;
486 e_x_loc[2] = u_2_parallel[2] * inv_u_2_parallel_norm;
487 }
488 T u2 = e_x_loc[0] * u_2[0] +
489 e_x_loc[1] * u_2[1] +
490 e_x_loc[2] * u_2[2]; //scalar Product of the basis vector e_x_loc and u_2
491
492 // STEP 5 : u1 can only be 0 if u2==0 for musker profile
493 //T* tau_w = new T [1];
494 T tau_w[1] = {0.};
495 T u1[1];
496 if (u2 != 0) {
497
498 T rho_phy = _converter.getPhysDensity(); // [kg/m³] SI units
499 T rho_lat = _converter.getLatticeDensity(rho_phy); // [kg/m³] DESCRIPTOR units
500 T Mol_Lat_Nu = _converter.getLatticeViscosity(); // [m²/s] DESCRIPTOR units
501
502 if (_wallFunctionParam.wallProfile == 0) {
503 // Musker profile for boundary node using the molecular viscosity and the characteristic density in lattice units
504 Musker<T,T> musker_u2_Lat_MolVis(Mol_Lat_Nu, y_2, rho_lat);
505 util::Newton1D<T> newton(musker_u2_Lat_MolVis,u2);
506 T tau_w_guess = cell.template getField<descriptors::TAU_W>();
507 tau_w[0] = newton.solve(tau_w_guess,false); // Computation of local wall shear stress
508
509 if (std::isnan(tau_w[0])) {
510 //OstreamManager clout(std::cout, "First NAN");
511 //clout << "Position = [" << x << ", " << y << "," << z << "]" << std::endl;
512 //clout << "Tau guess = " << tau_w_guess << std::endl;
513 tau_w[0] = newton.solve(0.,false);
514 }
515
516 if (std::isnan(tau_w[0])) {
517 OstreamManager clout(std::cout, "Second NAN");
518 clout << "Position = [" << x << ", " << y << "," << z << "]" << std::endl;
519 clout << "Tau guess = " << tau_w_guess << std::endl;
520 tau_w[0] = 0.;
521 }
522
523 // compute u1 - boundary node
524 Musker<T,T> musker_u1_Lat_MolVis(Mol_Lat_Nu, y_1, rho_lat);
525 musker_u1_Lat_MolVis(u1,tau_w);
526
527 }
528 else {
529 // Werner and Wengle profile for boundary node using the molecular viscosity and the characteristic density in lattice units
530 PowerLawProfile<T,T> PLP (Mol_Lat_Nu, u2, y_2, y_1, rho_lat);
531 T input [1];
532 T output[2];
533 PLP(output,input);
534 tau_w[0] = output[0];
535 u1[0] = output[1];
536 }
537 }
538 else {
539 u1[0] = 0.;
540 }
541 // save tau_w for next step
542 cell.template setField<descriptors::TAU_W>(&(tau_w[0]));
543 // STEP 6 : compute velocity vector at the boundary
544 u[0] = e_x_loc[0] * u1[0];
545 u[1] = e_x_loc[1] * u1[0];
546 u[2] = e_x_loc[2] * u1[0];
547
548 if (_wallFunctionParam.bodyForce) {
549 auto force = cell.template getFieldPointer<descriptors::FORCE>();
550 for (int iDim=0; iDim<DESCRIPTOR::d; ++iDim) {
551 u[iDim] -= force[iDim] / 2.;
552 }
553 }
554}
555
557template<typename T, typename DESCRIPTOR>
558void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeVanDriestTauEff(T y_bc, T tau_w, T u_bc, T u_1, T u_2, T& tau_eff)
559{
560
561 T rho_phy = _converter.getPhysDensity(); // [kg/m³] SI units
562 T rho = _converter.getLatticeDensity(rho_phy); // [kg/m³] DESCRIPTOR units
563 T nu_mol = _converter.getLatticeViscosity(); // [m²/s] DESCRIPTOR units
564
565 T y_plus = y_bc*util::sqrt(tau_w/rho)/nu_mol;
566 T uxdz_abs = util::abs(fd::boundaryGradient(u_bc, u_1, u_2));
567
568 T vanDriest = 1 - util::exp(-y_plus / 26.);
569
570 T nu_turb = util::pow(_wallFunctionParam.vonKarman * y_bc * vanDriest, 2.) * uxdz_abs;
571
572 T nu_eff = nu_turb + nu_mol;
573
574 tau_eff = (3. * nu_eff) + 0.5; // Rectangular integration rule - Collision
575}
577
578template<typename T, typename DESCRIPTOR>
579void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::ComputeTauEff(BlockLattice<T,DESCRIPTOR>& blockLattice, Cell<T,DESCRIPTOR>& cell, int x, int y, int z, T u_bc[DESCRIPTOR::d])
580{
581
582 T u_z[3];
583 T u_z2[3];
584 blockLattice.get(x - discreteNormalX, y - discreteNormalY, z - discreteNormalZ).computeU(u_z);
585 blockLattice.get(x - 2*discreteNormalX, y - 2*discreteNormalY, z - 2*discreteNormalZ).computeU(u_z2);
586
587 T y_bc = 0.5; // Default half-way Bounce-Back distance - [m] DESCRIPTOR units
588 if (_wallFunctionParam.latticeWalldistance > 0.) {
589 y_bc = _wallFunctionParam.latticeWalldistance; // [m] DESCRIPTOR units
590 }
591 T normal_norm = util::sqrt(discreteNormalX * discreteNormalX +
592 discreteNormalY * discreteNormalY +
593 discreteNormalZ * discreteNormalZ); // l2 norm : magnitude of the vector
594 y_bc *= normal_norm;
595 T tau_w = cell.template getField<descriptors::TAU_W>();
596
597 T tau_eff;
598
599 if (_wallFunctionParam.curved == true || orientation == 0) {
600 T inv_normal_norm = 1. / normal_norm;
601 T normal[3];
602 normal[0] = discreteNormalX * (inv_normal_norm);
603 normal[1] = discreteNormalY * (inv_normal_norm);
604 normal[2] = discreteNormalZ * (inv_normal_norm);
605
606 T u_2_parallel[3];
607 T sP_normal_u2 = normal[0] * u_z[0] +
608 normal[1] * u_z[1] +
609 normal[2] * u_z[2]; //scalar Product of the normal and u2
610
611 u_2_parallel[0] = u_z[0] - sP_normal_u2 * normal[0];
612 u_2_parallel[1] = u_z[1] - sP_normal_u2 * normal[1];
613 u_2_parallel[2] = u_z[2] - sP_normal_u2 * normal[2];
614
615 T u_2_parallel_norm = util::sqrt(u_2_parallel[0] * u_2_parallel[0] +
616 u_2_parallel[1] * u_2_parallel[1] +
617 u_2_parallel[2] * u_2_parallel[2]); // l2 norm : magnitude of the vector
618
619 T e_x_loc[3];
620 T inv_u_2_parallel_norm = (1. /u_2_parallel_norm );
621
622 if (u_2_parallel_norm != 0) {
623 e_x_loc[0] = u_2_parallel[0] * inv_u_2_parallel_norm;
624 e_x_loc[1] = u_2_parallel[1] * inv_u_2_parallel_norm;
625 e_x_loc[2] = u_2_parallel[2] * inv_u_2_parallel_norm;
626 }
627 else {
628 e_x_loc[0] = T(-discreteNormalX);
629 e_x_loc[1] = T(-discreteNormalY);
630 e_x_loc[2] = T(-discreteNormalZ);
631 }
632
633 T u1;
634 T u2;
635 T u3;
636 u1 = e_x_loc[0] * u_bc[0] +
637 e_x_loc[1] * u_bc[1] +
638 e_x_loc[2] * u_bc[2];
639 u2 = e_x_loc[0] * u_z[0] +
640 e_x_loc[1] * u_z[1] +
641 e_x_loc[2] * u_z[2];
642 u3 = e_x_loc[0] * u_z2[0] +
643 e_x_loc[1] * u_z2[1] +
644 e_x_loc[2] * u_z2[2];
645
646 computeVanDriestTauEff(y_bc, tau_w, u1, u2, u3, tau_eff);
647 }
648 else {
649 computeVanDriestTauEff(y_bc, tau_w, u_bc[direction], u_z[direction], u_z2[direction], tau_eff);
650 }
651 cell.template setField<descriptors::TAU_EFF>(&(tau_eff));
652
653}
654
655template<typename T, typename DESCRIPTOR>
656void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::ComputeRhoWall(BlockLattice<T,DESCRIPTOR>& blockLattice, Cell<T,DESCRIPTOR>& cell, int x, int y, int z, T u_bc[DESCRIPTOR::d], T& rho_bc)
657{
658
660 if (_wallFunctionParam.rhoMethod == 0) {
662 // The code have been copied from VelocityBM -> compute rho method. This implementation it is only valid for straight boundaries
663 T u_bc_Perpendicular = 0.;
664 if (_wallFunctionParam.bodyForce) {
665 T u_bc_tmp[DESCRIPTOR::d];
666 auto force = cell.template getFieldPointer<descriptors::FORCE>();
667 for (int iDim = 0; iDim<DESCRIPTOR::d; ++iDim) {
668 u_bc_tmp[iDim] = u_bc[iDim] - force[iDim] / 2.;
669 }
670 for (int iDim = 0; iDim<DESCRIPTOR::d; ++iDim) {
671 u_bc_Perpendicular += u_bc_tmp[iDim]*unit_normal[iDim];
672 }
673 }
674 else {
675 for (int iDim = 0; iDim<DESCRIPTOR::d; ++iDim) {
676 u_bc_Perpendicular += u_bc[iDim]*unit_normal[iDim];
677 }
678 }
679
680 T rhoOnWall = T();
681 for (unsigned fIndex=0; fIndex<onWallIndices.size(); ++fIndex) {
682 rhoOnWall += cell[onWallIndices[fIndex]];
683 }
684
685 T rhoNormal = T();
686 for (unsigned fIndex=0; fIndex<normalIndices.size(); ++fIndex) {
687 rhoNormal += cell[normalIndices[fIndex]];
688 }
689
690 rho_bc = ((T)2*rhoNormal+rhoOnWall+(T)1) /
691 ((T)1+u_bc_Perpendicular);
692
693 }
694 else if (_wallFunctionParam.rhoMethod == 1) { // Extrapolation Fluid
695
696 rho_bc = blockLattice.get(x - discreteNormalX, y - discreteNormalY, z - discreteNormalZ).computeRho();
697
698 }
699 else if (_wallFunctionParam.rhoMethod == 2) { // constant rho
700 rho_bc = 1.;
701 }
702}
703
705template<typename T, typename DESCRIPTOR>
707{
708 T sum_cell_fneq = 0.;
709 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
710 sum_cell_fneq += fneq_bc[iPop];
711 }
713 int iPi = 0;
714 for (int iAlpha=0; iAlpha < DESCRIPTOR::d; ++iAlpha) {
715 for (int iBeta=iAlpha; iBeta < DESCRIPTOR::d; ++iBeta) {
716 pi_bc[iPi] = T();
717 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
718 pi_bc[iPi] += descriptors::c<DESCRIPTOR>(iPop)[iAlpha]*descriptors::c<DESCRIPTOR>(iPop)[iBeta]*fneq_bc[iPop];
719 }
720 if (iAlpha==iBeta) {
721 pi_bc[iPi] -= (1./descriptors::invCs2<T,DESCRIPTOR>())*sum_cell_fneq;
722 }
723 ++iPi;
724 }
725 }
726
727 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
728 fneq_bc[iPop] = equilibrium<DESCRIPTOR>::template fromPiToFneq<T>(iPop, pi_bc);
729 }
730}
731
732template<typename T, typename DESCRIPTOR>
733void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeFneqRNEBB(Cell<T,DESCRIPTOR>& cell, T u_bc[DESCRIPTOR::d], T rho_bc, T fneq_bc[DESCRIPTOR::q])
734{
735 Dynamics<T,DESCRIPTOR>* dynamics = cell.getDynamics();
736 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
737 fneq_bc[iPop] = cell[iPop] - dynamics -> computeEquilibrium(iPop,rho_bc,u_bc);
738 }
739 for (unsigned fIndex=0; fIndex<normalInwardsIndices.size(); ++fIndex) {
740 fneq_bc[normalInwardsIndices[fIndex]] = fneq_bc[descriptors::opposite<DESCRIPTOR>(normalInwardsIndices[fIndex])];
741 }
742
743 computeRFneqfromFneq(fneq_bc);
744}
745
746template<typename T, typename DESCRIPTOR>
747void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeFneqENeq(BlockLattice<T,DESCRIPTOR>& blockLattice, Cell<T,DESCRIPTOR>& cell, int x, int y, int z, T u_bc[DESCRIPTOR::d], T rho_bc, T fneq_bc[DESCRIPTOR::q])
748{
749 Cell<T,DESCRIPTOR> cell_fluid = blockLattice.get(x - discreteNormalX, y - discreteNormalY, z - discreteNormalZ);
750 Dynamics<T,DESCRIPTOR>* dynamics_fluid = cell_fluid.getDynamics();
751 T rho_fluid, u_fluid[DESCRIPTOR::d];
752 cell_fluid.computeRhoU(rho_fluid,u_fluid);
753 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
754 fneq_bc[iPop] = cell_fluid[iPop] - dynamics_fluid -> computeEquilibrium(iPop,rho_fluid,u_fluid);
755 }
756}
757
758template<typename T, typename DESCRIPTOR>
759void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::computeFneqRSOFD(BlockLattice<T,DESCRIPTOR>& blockLattice, Cell<T,DESCRIPTOR>& cell, int x, int y, int z, T u_bc[DESCRIPTOR::d], T rho_bc, T fneq_bc[DESCRIPTOR::q])
760{
762 T Velocity_Grad[DESCRIPTOR::d][DESCRIPTOR::d];
763 computeVelocityGradient(blockLattice, x, y, z, u_bc, Velocity_Grad);
764 // Creation of strain Rate
765 T tau_eff = cell.template getField<descriptors::TAU_EFF>();
766 T Factor = -1.*tau_eff*rho_bc/descriptors::invCs2<T,DESCRIPTOR>();
767 int iPi = 0;
768 for (int Alpha=0; Alpha<DESCRIPTOR::d; ++Alpha) {
769 for (int Beta=Alpha; Beta<DESCRIPTOR::d; ++Beta) {
770 pi_bc[iPi] = Factor*(Velocity_Grad[Alpha][Beta] + Velocity_Grad[Beta][Alpha]);
771 ++iPi;
772 }
773 }
774
775 if (_wallFunctionParam.bodyForce) {
776 // Force tensor
777 //Creation of body force tensor (rho/2.)*(G_alpha*U_beta + U_alpha*G_Beta)
778 auto F = cell.template getField<descriptors::FORCE>();
779
780 iPi = 0;
781 for (int Alpha=0; Alpha<DESCRIPTOR::d; ++Alpha) {
782 for (int Beta=Alpha; Beta<DESCRIPTOR::d; ++Beta) {
783 pi_bc[iPi] += -1.*(rho_bc/2.)*(F[Alpha]*u_bc[Beta] + u_bc[Alpha]*F[Beta]);
784 ++iPi;
785 }
786 }
787 }
788
789 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
790 fneq_bc[iPop] = equilibrium<DESCRIPTOR>::template fromPiToFneq<T>(iPop, pi_bc);
791 }
792
793}
794
796
797template<typename T, typename DESCRIPTOR>
798void WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>::ComputeFneqWall(BlockLattice<T,DESCRIPTOR>& blockLattice, Cell<T,DESCRIPTOR>& cell, int x, int y, int z, T u_bc[DESCRIPTOR::d], T rho_bc, T fneq_bc[DESCRIPTOR::q])
799{
800
801 //regularized NEBB (Latt)
802 if (_wallFunctionParam.fneqMethod == 0) {
803 computeFneqRNEBB(cell, u_bc, rho_bc, fneq_bc);
804
805 //extrapolation NEQ (Guo Zhaoli)
806 }
807 else if (_wallFunctionParam.fneqMethod == 1) {
808 computeFneqENeq(blockLattice, cell, x, y, z, u_bc, rho_bc, fneq_bc);
809
810 //Regularized second order finite Differnce
811 }
812 else if (_wallFunctionParam.fneqMethod == 2) {
813 computeFneqRSOFD(blockLattice, cell, x, y, z, u_bc, rho_bc, fneq_bc);
814
815 //equilibrium scheme
816 }
817 else if (_wallFunctionParam.fneqMethod == 3) {
818 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
819 fneq_bc[iPop] = 0.;
820 }
821 }
822
823}
824
825template<typename T, typename DESCRIPTOR>
827{
828 Cell<T,DESCRIPTOR> cell_bc = blockLattice.get(x,y,z);
829 T rho_bc = 0.;
830 T u_bc[DESCRIPTOR::d];
831 T fneq_bc[DESCRIPTOR::q];
832
833 // Computation of boundary velocity from the wall function
834 ComputeUWall(blockLattice, x, y, z, u_bc);
835
836 if (_wallFunctionParam.useVanDriest) {
837 // Computation of effective relaxation time from the vanDriest damping function
838 ComputeTauEff(blockLattice, cell_bc, x, y, z, u_bc);
839 }
840
841 // Computation of density at the boundary node
842 ComputeRhoWall(blockLattice, cell_bc, x, y, z, u_bc, rho_bc);
843
844 // Computation of the second-order moment of non-equilibrium distribution function
845 ComputeFneqWall(blockLattice, cell_bc, x, y, z, u_bc, rho_bc, fneq_bc);
846
847 // Computation of the particle distribution functions according to the regularized formula
848 Dynamics<T,DESCRIPTOR>* dynamics_bc = cell_bc.getDynamics();
849
850 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
851 cell_bc[iPop] = dynamics_bc -> computeEquilibrium(iPop,rho_bc,u_bc) + fneq_bc[iPop];
852 if (std::isnan(cell_bc[iPop])) {
853 OstreamManager clout(std::cout, "Slip Musker Profile");
854 clout << "Musker Profile Computation" << std::endl;
855 clout << "Position = [" << x << ", " << y << "," << z << "]" << std::endl;
856 clout << "Normal outwards = [" << discreteNormalX << ", " << discreteNormalY << "," << discreteNormalZ << "]" << std::endl;
857 clout << "Velocity at boundary u_bc = [" << u_bc[0] << ", " << u_bc[1] << "," << u_bc[2] << "]" << std::endl;
858 clout << "Density at boundary rho_bc = " << rho_bc << std::endl;
859 exit(1);
860 }
861 }
862
863}
864
865template<typename T, typename DESCRIPTOR>
867{
868 processSubDomain(blockLattice, x0, x1, y0, y1, z0, z1);
869}
870
872
873template<typename T, typename DESCRIPTOR>
875 std::vector<int> discreteNormal, std::vector<int> missingIndices,
877 IndicatorF3D<T>* geoIndicator)
878 : PostProcessorGenerator3D<T,DESCRIPTOR>(x0, x1, y0, y1, z0, z1),
879 _blockGeometryStructure(blockGeometryStructure),
880 _discreteNormal(discreteNormal), _missingIndices(missingIndices),
881 _converter(converter), _wallFunctionParam(wallFunctionParam), _geoIndicator(geoIndicator)
882{ }
883
884template<typename T, typename DESCRIPTOR>
887{
888 return new WallFunctionBoundaryProcessor3D<T,DESCRIPTOR>( this->x0, this->x1, this->y0, this->y1, this->z0, this->z1,
889 _blockGeometryStructure, _discreteNormal, _missingIndices,
890 _converter, _wallFunctionParam, _geoIndicator);
891}
892
893template<typename T, typename DESCRIPTOR>
896{
897 return new WallFunctionBoundaryProcessorGenerator3D<T,DESCRIPTOR> (this->x0, this->x1, this->y0, this->y1, this->z0, this->z1,
898 _blockGeometryStructure, _discreteNormal, _missingIndices,
899 _converter, _wallFunctionParam, _geoIndicator);
900}
901
902} // namespace olb
903
904#endif
AnalyticalF are applications from DD to XD, where X is set by the constructor.
Representation of a block geometry.
Platform-abstracted block lattice for external access and inter-block interaction.
Cell< T, DESCRIPTOR > get(CellID iCell)
Get Cell interface for index iCell.
Highest-level interface to Cell data.
Definition cell.h:148
Dynamics< T, DESCRIPTOR > * getDynamics()
Get a pointer to the dynamics.
Definition cell.hh:174
void computeRhoU(T &rho, T u[descriptors::d< DESCRIPTOR >()]) const
Compute fluid velocity and particle density on the cell.
Definition cell.hh:232
std::string & getName()
read and write access to name
Definition genericF.hh:51
IndicatorF3D is an application from .
bool operator()(T output[], const S tau_w[])
has to be implemented for 'every' derived class
class for marking output with some text
std::string & getName()
read and write access to name
PowerLawProfile(T nu, T u2, T y2, T y1, T rho)
bool operator()(T output[], const S tau_w[])
has to be implemented for 'every' derived class
Conversion between physical and lattice units, as well as discretization.
WallFunctionBoundaryProcessor3D(int x0, int x1, int y0, int y1, int z0, int z1, BlockGeometry< T, 3 > &blockGeometryStructure, std::vector< int > discreteNormal, std::vector< int > missingIndices, UnitConverter< T, DESCRIPTOR > const &converter, wallFunctionParam< T > const &wallFunctionParam, IndicatorF3D< T > *geoIndicator)
virtual void processSubDomain(BlockLattice< T, DESCRIPTOR > &blockLattice, int x0_, int x1_, int y0_, int y1_, int z0_, int z1_)
Execute post-processing step on a sublattice.
virtual void ComputeWallFunction(BlockLattice< T, DESCRIPTOR > &blockLattice, int x, int y, int z)
virtual void process(BlockLattice< T, DESCRIPTOR > &blockLattice)
Execute post-processing step.
PostProcessorGenerator3D< T, DESCRIPTOR > * clone() const override
WallFunctionBoundaryProcessorGenerator3D(int x0, int x1, int y0, int y1, int z0, int z1, BlockGeometry< T, 3 > &blockGeometryStructure, std::vector< int > discreteNormal, std::vector< int > missingIndices, UnitConverter< T, DESCRIPTOR > const &converter, wallFunctionParam< T > const &wallFunctionParam, IndicatorF3D< T > *geoIndicator)
PostProcessor3D< T, DESCRIPTOR > * generate() const override
1D Newton simple scheme
constexpr T BSGradient(T u_0, T u_1, T u_2)
Backward second-order first derivative.
constexpr T centralGradient(T u_p1, T u_m1) any_platform
Second-order central gradient (u_p1 = u(x+1))
constexpr T FSGradient(T u_0, T u_1, T u_2)
Forward second-order first derivative.
constexpr T boundaryGradient(T u_0, T u_1, T u_2) any_platform
Second-order asymmetric gradient (u_1 = u(x+1))
cpu::simd::Pack< T > sqrt(cpu::simd::Pack< T > value)
Definition pack.h:100
ADf< T, DIM > abs(const ADf< T, DIM > &a)
Definition aDiff.h:1019
bool distance(S &distance, const Vector< S, D > &origin, const Vector< S, D > &direction, S precision, S pitch, F1 isInside, F2 isInsideBoundingBox)
ADf< T, DIM > log(const ADf< T, DIM > &a)
Definition aDiff.h:475
cpu::simd::Pack< T > pow(cpu::simd::Pack< T > base, cpu::simd::Pack< T > exp)
Definition pack.h:112
bool intersect(int x0, int x1, int y0, int y1, int x0_, int x1_, int y0_, int y1_, int &newX0, int &newX1, int &newY0, int &newY1)
Definition util.h:89
ADf< T, DIM > atan(const ADf< T, DIM > &a)
Definition aDiff.h:614
ADf< T, DIM > exp(const ADf< T, DIM > &a)
Definition aDiff.h:455
Top level namespace for all of OpenLB.
constexpr T norm(const ScalarVector< T, D, IMPL > &a)
Euclidean vector norm.
Interface for per-cell dynamics.
Definition interface.h:54
Compute number of elements of a symmetric d-dimensional tensor.
Definition util.h:210
Set of functions commonly used in LB computations – header file.