ParserUBB.h 23.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
//======================================================================================================================
//
//  This file is part of waLBerla. waLBerla is free software: you can 
//  redistribute it and/or modify it under the terms of the GNU General Public
//  License as published by the Free Software Foundation, either version 3 of 
//  the License, or (at your option) any later version.
//  
//  waLBerla is distributed in the hope that it will be useful, but WITHOUT 
//  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 
//  FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License 
//  for more details.
//  
//  You should have received a copy of the GNU General Public License along
//  with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file ParserUBB.h
//! \ingroup lbm
//! \author Michael Kuron <mkuron@icp.uni-stuttgart.de>
//
//======================================================================================================================

#pragma once

#include "lbm/field/DensityAndVelocity.h"
#include "lbm/field/PdfField.h"
#include "lbm/refinement/TimeTracker.h"

#include "boundary/Boundary.h"

#include "core/DataTypes.h"
#include "core/cell/CellInterval.h"
#include "core/config/Config.h"
#include "core/debug/Debug.h"
#include "core/math/Vector3.h"
#include "core/math/ParserOMP.h"

#include "field/FlagUID.h"

#include "stencil/Directions.h"

#include <vector>



namespace walberla {
namespace lbm {



50
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce = false, bool StoreForce = false >
51
52
53
54
55
class ParserUBB : public Boundary<flag_t>
{
   typedef lbm::PdfField< LatticeModel_T >   PDFField;
   typedef typename LatticeModel_T::Stencil  Stencil;

56
57
   typedef GhostLayerField< Vector3<real_t>, 1 > ForceField;

58
59
60
61
62
63
64
65
public:

   static const bool threadsafe = true;

   class Parser : public BoundaryConfiguration
   {
   public:
      inline Parser( const Config::BlockHandle & config );
Michael Kuron's avatar
Michael Kuron committed
66
      inline Parser( const std::array< std::string, 3 > & equations );
67
68
      Vector3< real_t > operator()( const Vector3< real_t > & x, const real_t t ) const;
      Vector3< real_t > operator()( const Vector3< real_t > & x ) const;
69
      bool isTimeDependent() const { return timeDependent_; }
70
      const std::array< std::string, 3 > & equations() const { return equations_; }
71
72
73

   private:
      std::array< math::FunctionParserOMP, 3 > parsers_;
74
      std::array< std::string, 3 > equations_;
75
      bool timeDependent_;
76
77
   }; // class Parser

78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
   // constant velocity class is for API compatibility with UBB
   class Velocity : public BoundaryConfiguration {
   public:
             Velocity( const Vector3< real_t > & _velocity ) : velocity_( _velocity ) {}
             Velocity( const real_t _x, const real_t _y, const real_t _z ) : velocity_(_x,_y,_z) {}
      inline Velocity( const Config::BlockHandle & config );

      const Vector3< real_t > & velocity() const { return velocity_; }

      const real_t & x() const { return velocity_[0]; }
      const real_t & y() const { return velocity_[1]; }
      const real_t & z() const { return velocity_[2]; }

      Vector3< real_t > & velocity() { return velocity_; }

      real_t & x() { return velocity_[0]; }
      real_t & y() { return velocity_[1]; }
      real_t & z() { return velocity_[2]; }

   private:
      Vector3< real_t > velocity_;
   }; // class Velocity

101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
   typedef GhostLayerField< shared_ptr<Parser>, 1 > ParserField;
   typedef GhostLayerField< Vector3<real_t>, 1 >    VelocityField;

   static shared_ptr<Parser> createConfiguration( const Config::BlockHandle & config )
      { return make_shared<Parser>( config ); }



   ParserUBB( const BoundaryUID & boundaryUID, const FlagUID & uid, PDFField * const pdfField,
              FlagField<flag_t> * const flagField, const shared_ptr< TimeTracker > & timeTracker,
              const uint_t level, const AABB & aabb );
   ParserUBB( const BoundaryUID & boundaryUID, const FlagUID & uid, PDFField * const pdfField,
              FlagField<flag_t> * const flagField,  const uint_t level, const AABB & aabb );

   shared_ptr< TimeTracker > getTimeTracker() { return timeTracker_; }

   void pushFlags( std::vector< FlagUID > & uids ) const { uids.push_back( uid_ ); }

   void beforeBoundaryTreatment()
   {
      if( timeTracker_ )
         time_ = timeTracker_->getTime( level_ );
123
124
125

      if (StoreForce)
         force_->setWithGhostLayer( Vector3<real_t>() );
126
127
128
129
   }
   void  afterBoundaryTreatment() const {}

   template< typename Buffer_T >
130
   inline void packCell( Buffer_T &, const cell_idx_t, const cell_idx_t, const cell_idx_t ) const;
131
132

   template< typename Buffer_T >
133
   inline void registerCell( Buffer_T &, const flag_t, const cell_idx_t, const cell_idx_t, const cell_idx_t );
134
135
136
137
138
139

   inline void registerCell( const flag_t, const cell_idx_t x, const cell_idx_t y, const cell_idx_t z, const BoundaryConfiguration & parser );
   inline void registerCells( const flag_t, const CellInterval & cells, const BoundaryConfiguration & parser );
   template< typename CellIterator >
   inline void registerCells( const flag_t, const CellIterator & begin, const CellIterator & end, const BoundaryConfiguration & parser );

140
   inline void unregisterCell( const flag_t, const cell_idx_t, const cell_idx_t, const cell_idx_t ) const;
141
142
143
144
145
146
147
148
149

#ifndef NDEBUG
   inline void treatDirection( const cell_idx_t  x, const cell_idx_t  y, const cell_idx_t  z, const stencil::Direction dir,
                               const cell_idx_t nx, const cell_idx_t ny, const cell_idx_t nz, const flag_t mask );
#else
   inline void treatDirection( const cell_idx_t  x, const cell_idx_t  y, const cell_idx_t  z, const stencil::Direction dir,
                               const cell_idx_t nx, const cell_idx_t ny, const cell_idx_t nz, const flag_t /*mask*/ );
#endif

150
151
152
   inline Vector3<real_t> getValue( const cell_idx_t x, cell_idx_t y, cell_idx_t z ) const;
   inline Vector3<real_t> getValue( const cell_idx_t x, cell_idx_t y, cell_idx_t z, real_t t ) const;

153
154
155
156
157
158
   const typename ForceField::value_type & getForce( const cell_idx_t x, cell_idx_t y, cell_idx_t z ) const
   {
      static_assert(StoreForce, "this member function is only available if the fourth template argument on the class is true");
      return force_->get(x,y,z);
   }

159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
private:

   const FlagUID uid_;

   PDFField * const pdfField_;

   Vector3< real_t > origin_;
   Vector3< real_t > dx_;

   // required to keep track of the simulation time
   shared_ptr< TimeTracker > timeTracker_; // -> addPostBoundaryHandlingVoidFunction (when used with refinement time step)
   real_t time_;
   uint_t level_;

   shared_ptr<ParserField> parserField_;
   shared_ptr<VelocityField> velocityField_;
175
   shared_ptr<ForceField> force_;
176
177
178

}; // class ParserUBB

179
180
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
inline ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::Parser::Parser( const Config::BlockHandle & config )
181
: parsers_(), equations_(), timeDependent_( false )
182
183
184
185
186
187
{
   if( !config )
      return;

   if( config.isDefined( "x" ) )
   {
188
189
      equations_[0] = config.getParameter<std::string>( "x" );
      parsers_[0].parse( equations_[0] );
190
      if( parsers_[0].symbolExists( "t" ) )
191
         timeDependent_ = true;
192
193
194
   }
   if( config.isDefined( "y" ) )
   {
195
196
      equations_[1] = config.getParameter<std::string>( "y" );
      parsers_[1].parse( equations_[1] );
197
      if( parsers_[1].symbolExists( "t" ) )
198
         timeDependent_ = true;
199
200
201
   }
   if( config.isDefined( "z" ) )
   {
202
203
204
205
206
207
208
      equations_[2] = config.getParameter<std::string>( "z" );
      parsers_[2].parse( equations_[2] );
      if( parsers_[2].symbolExists( "t" ) )
         timeDependent_ = true;
   }
}

209
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
Michael Kuron's avatar
Michael Kuron committed
210
inline ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::Parser::Parser( const std::array< std::string, 3 > & equations )
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
: parsers_(), equations_( equations ), timeDependent_( false )
{
   if( equations_[0].length() > 0 )
   {
      parsers_[0].parse( equations_[0] );
      if( parsers_[0].symbolExists( "t" ) )
         timeDependent_ = true;
   }
   if( equations_[1].length() > 0 )
   {
      parsers_[1].parse( equations_[1] );
      if( parsers_[1].symbolExists( "t" ) )
         timeDependent_ = true;
   }
   if( equations_[2].length() > 0 )
   {
      parsers_[2].parse( equations_[2] );
228
      if( parsers_[2].symbolExists( "t" ) )
229
         timeDependent_ = true;
230
231
232
   }
}

233
234
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::Velocity::Velocity( const Config::BlockHandle & config  )
235
236
237
238
239
240
{
   velocity_[0] = ( config && config.isDefined( "x" ) ) ? config.getParameter<real_t>( "x" ) : real_c(0.0);
   velocity_[1] = ( config && config.isDefined( "y" ) ) ? config.getParameter<real_t>( "y" ) : real_c(0.0);
   velocity_[2] = ( config && config.isDefined( "z" ) ) ? config.getParameter<real_t>( "z" ) : real_c(0.0);
}

241
242
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
Vector3< real_t > ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::Parser::operator()( const Vector3< real_t > & x, const real_t t ) const
243
244
245
246
247
248
249
250
251
252
253
254
255
256
{
   std::map< std::string, double > symbols;
   symbols["x"] = x[0];
   symbols["y"] = x[1];
   symbols["z"] = x[2];
   symbols["t"] = t;

   Vector3< real_t > v;
   v[0] = parsers_[0].evaluate( symbols );
   v[1] = parsers_[1].evaluate( symbols );
   v[2] = parsers_[2].evaluate( symbols );
   return v;
}

257
258
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
Vector3< real_t > ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::Parser::operator()( const Vector3< real_t > & x ) const
259
{
260
   WALBERLA_ASSERT( !timeDependent_ );
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275

   std::map< std::string, double > symbols;
   symbols["x"] = x[0];
   symbols["y"] = x[1];
   symbols["z"] = x[2];

   Vector3< real_t > v;
   v[0] = parsers_[0].evaluate( symbols );
   v[1] = parsers_[1].evaluate( symbols );
   v[2] = parsers_[2].evaluate( symbols );
   return v;
}



276
277
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
inline ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::ParserUBB( const BoundaryUID & boundaryUID, const FlagUID & uid, PDFField * const pdfField,
278
279
                                                                                   FlagField<flag_t> * const flagField, const shared_ptr< TimeTracker > & timeTracker,
                                                                                   const uint_t level, const AABB & aabb )
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
   : Boundary<flag_t>( boundaryUID ), uid_( uid ), pdfField_( pdfField ), timeTracker_( timeTracker ), time_( real_t(0) ), level_( level )
{
   WALBERLA_ASSERT_NOT_NULLPTR( pdfField_ );
   dx_[0] = aabb.xSize() / real_c( pdfField_->xSize() );
   dx_[1] = aabb.ySize() / real_c( pdfField_->ySize() );
   dx_[2] = aabb.zSize() / real_c( pdfField_->zSize() );
   origin_[0] = aabb.xMin() + real_c(0.5) * dx_[0];
   origin_[1] = aabb.yMin() + real_c(0.5) * dx_[1];
   origin_[2] = aabb.zMin() + real_c(0.5) * dx_[2];

   if(flagField != NULL)
   {
      parserField_   = make_shared<ParserField>  ( pdfField->xSize(), pdfField->ySize(), pdfField->zSize(), flagField->nrOfGhostLayers(), field::zyxf );
      velocityField_ = make_shared<VelocityField>( pdfField->xSize(), pdfField->ySize(), pdfField->zSize(), flagField->nrOfGhostLayers(), field::zyxf );
   }
   else
   {
      parserField_   = make_shared<ParserField>  ( pdfField->xSize(), pdfField->ySize(), pdfField->zSize(), pdfField->nrOfGhostLayers(),  field::zyxf );
      velocityField_ = make_shared<VelocityField>( pdfField->xSize(), pdfField->ySize(), pdfField->zSize(), pdfField->nrOfGhostLayers(),  field::zyxf );
   }
300
301
302

   if (StoreForce)
      force_ = make_shared<ForceField>( pdfField_->xSize(), pdfField_->ySize(), pdfField_->zSize(), pdfField_->nrOfGhostLayers(), field::zyxf );
303
304
305
306
}



307
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
308
template< typename Buffer_T >
309
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::packCell( Buffer_T & buffer, const cell_idx_t x, const cell_idx_t y, const cell_idx_t z ) const
310
311
312
313
314
315
316
317
318
319
320
321
322
323
{
   if( parserField_->get( x, y, z ) )
   {
      auto & eqs = parserField_->get( x, y, z )->equations();
      buffer << true << eqs[0] << eqs[1] << eqs[2];
   }
   else
   {
      buffer << false << velocityField_->get( x, y, z );
   }
}



324
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
325
template< typename Buffer_T >
326
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::registerCell( Buffer_T & buffer, const flag_t, const cell_idx_t x, const cell_idx_t y, const cell_idx_t z )
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
{
   bool isparser;
   buffer >> isparser;
   if( isparser )
   {
      std::array< std::string, 3> eqs;
      buffer >> eqs[0] >> eqs[1] >> eqs[2];

      auto p = make_shared<Parser>(eqs);
      parserField_->get( x, y, z ) = p;
   }
   else
   {
      buffer >> velocityField_->get( x, y, z );
      parserField_->get( x, y, z ) = nullptr;
   }
}



347
348
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
inline ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::ParserUBB( const BoundaryUID & boundaryUID, const FlagUID & uid, PDFField * const pdfField,
349
                                                                                   FlagField<flag_t> * const flagField, const uint_t level, const AABB & aabb )
350
351
352
353
354
   : ParserUBB( boundaryUID, uid, pdfField, flagField, nullptr, level, aabb )
{}



355
356
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::registerCell( const flag_t, const cell_idx_t x, const cell_idx_t y, const cell_idx_t z,
357
358
359
360
361
                                                                                             const BoundaryConfiguration & parser )
{
   WALBERLA_ASSERT_EQUAL( dynamic_cast< const Parser * >( &parser ), &parser );
   WALBERLA_ASSERT_NOT_NULLPTR( parserField_ );

362
363
364
365
366
367
368
   if( auto v = dynamic_cast< const Velocity * >( &parser ) )
   {
      velocityField_->get( x, y, z ) = v->velocity();
      parserField_->get( x, y, z ) = nullptr;
      return;
   }

369
370
371
372
   auto & p = dynamic_cast< const Parser & >( parser );

   if( p.isTimeDependent() )
   {
Michael Kuron's avatar
Michael Kuron committed
373
      parserField_->get( x, y, z ) = make_shared<Parser>( p.equations() );
374
375
376
377
378
379
380
381
382
383
384
385
386
   }
   else
   {
      const Vector3< real_t > pos( origin_[0] + real_c(x) * dx_[0],
                                   origin_[1] + real_c(y) * dx_[1],
                                   origin_[2] + real_c(z) * dx_[2] );
      velocityField_->get( x, y, z ) = p(pos);
      parserField_->get( x, y, z ) = nullptr;
   }
}



387
388
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::registerCells( const flag_t, const CellInterval & cells, const BoundaryConfiguration & parser )
389
390
391
392
{
   WALBERLA_ASSERT_EQUAL( dynamic_cast< const Parser * >( &parser ), &parser );
   WALBERLA_ASSERT_NOT_NULLPTR( parserField_ );

393
394
395
396
397
398
399
400
401
402
   if( auto v = dynamic_cast< const Velocity * >( &parser ) )
   {
      for( auto cell = parserField_->beginSliceXYZ( cells ); cell != parserField_->end(); ++cell )
      {
         velocityField_->get( cell.x(), cell.y(), cell.z() ) = v->velocity();
         *cell = nullptr;
      }
      return;
   }

403
404
405
406
   auto & p = dynamic_cast< const Parser & >( parser );

   if( p.isTimeDependent() )
   {
Michael Kuron's avatar
Michael Kuron committed
407
      auto shared_p = make_shared<Parser>( p.equations() );
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
      for( auto cell = parserField_->beginSliceXYZ( cells ); cell != parserField_->end(); ++cell )
         *cell = shared_p;
   }
   else
   {
      for( auto cell = parserField_->beginSliceXYZ( cells ); cell != parserField_->end(); ++cell )
      {
         const Vector3< real_t > pos( origin_[0] + real_c(cell.x()) * dx_[0],
                                      origin_[1] + real_c(cell.y()) * dx_[1],
                                      origin_[2] + real_c(cell.z()) * dx_[2] );
         velocityField_->get( cell.x(), cell.y(), cell.z() ) = p(pos);
         *cell = nullptr;
      }
   }
}



426
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
427
template< typename CellIterator >
428
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::registerCells( const flag_t, const CellIterator & begin, const CellIterator & end,
429
430
431
432
433
                                                                                              const BoundaryConfiguration & parser )
{
   WALBERLA_ASSERT_EQUAL( dynamic_cast< const Parser * >( &parser ), &parser );
   WALBERLA_ASSERT_NOT_NULLPTR( parserField_ );

434
435
436
437
438
439
440
441
442
443
   if( auto v = dynamic_cast< const Velocity * >( &parser ) )
   {
      for( auto cell = begin; cell != end; ++cell )
      {
         velocityField_->get( cell.x(), cell.y(), cell.z() ) = v->velocity();
         parserField_->get( cell->x(), cell->y(), cell->z() ) = nullptr;
      }
      return;
   }

444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
   auto & p = dynamic_cast< const Parser & >( parser );

   if( p.isTimeDependent() )
   {
      auto shared_p = make_shared<Parser>(p);
      for( auto cell = begin; cell != end; ++cell )
      {
         parserField_->get( cell->x(), cell->y(), cell->z() ) = shared_p;
      }
   }
   else
   {
      for( auto cell = begin; cell != end; ++cell )
      {
         const Vector3< real_t > pos( origin_[0] + real_c(cell->x()) * dx_[0],
                                      origin_[1] + real_c(cell->y()) * dx_[1],
                                      origin_[2] + real_c(cell->z()) * dx_[2] );
         velocityField_->get( cell->x(), cell->y(), cell->z() ) = p(pos);
         parserField_->get( cell->x(), cell->y(), cell->z() ) = nullptr;
      }
   }
}



469
470
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline void ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::unregisterCell( const flag_t, const cell_idx_t x, const cell_idx_t y, const cell_idx_t z ) const
471
472
473
474
475
476
{
   parserField_->get(x,y,z) = nullptr;
}



477
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce>
478
#ifndef NDEBUG
479
inline void ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::ParserUBB::treatDirection( const cell_idx_t  x, const cell_idx_t  y, const cell_idx_t  z, const stencil::Direction dir,
480
                                                                                                        const cell_idx_t nx, const cell_idx_t ny, const cell_idx_t nz, const flag_t mask )
481
#else
482
inline void ParserUBB<LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce>::ParserUBB::treatDirection( const cell_idx_t  x, const cell_idx_t  y, const cell_idx_t  z, const stencil::Direction dir,
483
                                                                                                        const cell_idx_t nx, const cell_idx_t ny, const cell_idx_t nz, const flag_t /*mask*/ )
484
485
486
487
488
489
490
491
492
#endif
   {
      WALBERLA_ASSERT_EQUAL( nx, x + cell_idx_c( stencil::cx[ dir ] ) );
      WALBERLA_ASSERT_EQUAL( ny, y + cell_idx_c( stencil::cy[ dir ] ) );
      WALBERLA_ASSERT_EQUAL( nz, z + cell_idx_c( stencil::cz[ dir ] ) );
      WALBERLA_ASSERT_UNEQUAL( mask & this->mask_, numeric_cast<flag_t>(0) );
      WALBERLA_ASSERT_EQUAL( mask & this->mask_, this->mask_ ); // only true if "this->mask_" only contains one single flag, which is the case for the
                                                                // current implementation of this boundary condition (ParserUBB)

493
494
      const real_t pdf_old = pdfField_->get( x, y, z, Stencil::idx[dir] );

495
496
497
      Vector3<real_t> velocity;
      if( parserField_->get(nx,ny,nz) )
      {
498
         WALBERLA_ASSERT_NOT_NULLPTR( getTimeTracker(), "A TimeTracker is needed for time-dependent equations" );
499
500
501
         const Vector3< real_t > pos( origin_[0] + real_c(nx) * dx_[0],
                                      origin_[1] + real_c(ny) * dx_[1],
                                      origin_[2] + real_c(nz) * dx_[2] );
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
         velocity = (*parserField_->get(nx,ny,nz))( pos, time_ );
      }
      else
      {
         velocity = velocityField_->get(nx,ny,nz);
      }

      if( LatticeModel_T::compressible )
      {
         const auto density  = pdfField_->getDensity(x,y,z);
         const auto vel = AdaptVelocityToExternalForce ? lbm::internal::AdaptVelocityToForce<LatticeModel_T>::get( x, y, z, pdfField_->latticeModel(), velocity, density ) :
                                                         velocity;

         pdfField_->get( nx, ny, nz, Stencil::invDirIdx(dir) ) = pdfField_->get( x, y, z, Stencil::idx[dir] ) -
                                                                 ( real_c(6) * density * real_c(LatticeModel_T::w[ Stencil::idx[dir] ]) *
                                                                    ( real_c(stencil::cx[ dir ]) * vel[0] +
                                                                      real_c(stencil::cy[ dir ]) * vel[1] +
                                                                      real_c(stencil::cz[ dir ]) * vel[2] ) );
      }
      else
      {
         const auto vel = AdaptVelocityToExternalForce ? lbm::internal::AdaptVelocityToForce<LatticeModel_T>::get( x, y, z, pdfField_->latticeModel(), velocity, real_t(1) ) :
                                                         velocity;

         pdfField_->get( nx, ny, nz, Stencil::invDirIdx(dir) ) = pdfField_->get( x, y, z, Stencil::idx[dir] ) -
                                                                 ( real_c(6) * real_c(LatticeModel_T::w[ Stencil::idx[dir] ]) *
                                                                    ( real_c(stencil::cx[ dir ]) * vel[0] +
                                                                      real_c(stencil::cy[ dir ]) * vel[1] +
                                                                      real_c(stencil::cz[ dir ]) * vel[2] ) );
      }
532
533
534
535
536
537
538
539
540

      if (StoreForce && pdfField_->isInInnerPart( Cell(x,y,z) ))
      {
         const real_t forceMEM = pdf_old + pdfField_->get( nx, ny, nz, Stencil::invDirIdx(dir) );
         Vector3<real_t> force( real_c( stencil::cx[dir] ) * forceMEM,
                                real_c( stencil::cy[dir] ) * forceMEM,
                                real_c( stencil::cz[dir] ) * forceMEM );
         force_->get( nx, ny, nz ) += force;
      }
541
542
543
   }


544

545
546
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline Vector3<real_t> ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::getValue( const cell_idx_t x, cell_idx_t y, cell_idx_t z ) const
547
548
549
550
{
   return getValue( x, y, z, time_ );
}

551
552
template< typename LatticeModel_T, typename flag_t, bool AdaptVelocityToExternalForce, bool StoreForce >
inline Vector3<real_t> ParserUBB< LatticeModel_T, flag_t, AdaptVelocityToExternalForce, StoreForce >::getValue( const cell_idx_t x, cell_idx_t y, cell_idx_t z, real_t t ) const
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
{
   Vector3<real_t> velocity;
   if( parserField_->get(x,y,z) )
   {
      WALBERLA_ASSERT_NOT_NULLPTR( getTimeTracker(), "A TimeTracker is needed for time-dependent equations" );
      const Vector3< real_t > pos( origin_[0] + real_c(x) * dx_[0],
                                   origin_[1] + real_c(y) * dx_[1],
                                   origin_[2] + real_c(z) * dx_[2] );
      velocity = (*parserField_->get(x,y,z))( pos, t );
   }
   else
   {
      velocity = velocityField_->get(x,y,z);
   }
   return velocity;
}


571
572
} // namespace lbm
} // namespace walberla