BioDynaMo  v1.05.119-a4ff3934
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
bdm::DiffusionGrid Class Referenceabstract

#include <diffusion_grid.h>

Inheritance diagram for bdm::DiffusionGrid:
[legend]
Collaboration diagram for bdm::DiffusionGrid:
[legend]

Public Member Functions

 DiffusionGrid ()=default
 
 DiffusionGrid (const TRootIOCtor *)
 
 DiffusionGrid (int substance_id, const std::string &substance_name, real_t dc, real_t mu, int resolution=10)
 
 ~DiffusionGrid () override=default
 
 DiffusionGrid (const DiffusionGrid &)=delete
 
DiffusionGridoperator= (const DiffusionGrid &)=delete
 
 DiffusionGrid (DiffusionGrid &&)=delete
 
DiffusionGridoperator= (DiffusionGrid &&)=delete
 
void Initialize () override
 
void Update () override
 
void Step (real_t dt) override
 
void Diffuse (real_t dt)
 
virtual void DiffuseWithClosedEdge (real_t dt)=0
 
virtual void DiffuseWithOpenEdge (real_t dt)=0
 
virtual void DiffuseWithDirichlet (real_t dt)=0
 
virtual void DiffuseWithNeumann (real_t dt)=0
 
virtual void DiffuseWithPeriodic (real_t dt)=0
 
void CalculateGradient ()
 
void RunInitializers ()
 
void ChangeConcentrationBy (const Real3 &position, real_t amount, InteractionMode mode=InteractionMode::kAdditive, bool scale_with_resolution=false)
 
void ChangeConcentrationBy (size_t idx, real_t amount, InteractionMode mode=InteractionMode::kAdditive, bool scale_with_resolution=false)
 See ChangeConcentrationBy(const Real3&, real_t, InteractionMode, bool). More...
 
real_t GetValue (const Real3 &position) const override
 Get the value of the scalar field at specified position. More...
 
real_t GetConcentration (const Real3 &position) const
 
real_t GetConcentration (const size_t idx) const
 Get the concentration at specified index. More...
 
Real3 GetGradient (const Real3 &position) const override
 
virtual void GetGradient (const Real3 &position, Real3 *gradient, bool normalize=true) const
 
std::array< uint32_t, 3 > GetBoxCoordinates (const Real3 &position) const
 Get the coordinates of the box at the specified position. More...
 
std::array< uint32_t, 3 > GetBoxCoordinates (const size_t idx) const
 
size_t GetBoxIndex (const std::array< uint32_t, 3 > &box_coord) const
 Calculates the box index at specified box coordinates. More...
 
size_t GetBoxIndex (const Real3 &position) const
 Calculates the box index of the substance at specified position. More...
 
std::array< size_t, 6 > GetNeighboringBoxes (size_t index) const
 
std::array< size_t, 6 > GetNeighboringBoxes (size_t index, const std::array< uint32_t, 3 > &box_coord) const
 
void SetDecayConstant (real_t mu)
 
real_t GetLastTimestep () const
 Return the last timestep dt that was used to run Diffuse(dt) More...
 
void SetUpperThreshold (real_t t)
 
real_t GetUpperThreshold () const
 
void SetLowerThreshold (real_t t)
 
real_t GetLowerThreshold () const
 
const real_tGetAllConcentrations () const
 
const real_tGetAllGradients () const
 
std::array< size_t, 3 > GetNumBoxesArray () const
 
size_t GetNumBoxes () const
 
real_t GetBoxLength () const
 
int GetSubstanceId () const
 
const std::string & GetSubstanceName () const
 
real_t GetDecayConstant () const
 
const int32_t * GetDimensionsPtr () const
 
std::array< int32_t, 6 > GetDimensions () const
 
std::array< int32_t, 3 > GetGridSize () const
 
const std::array< real_t, 7 > & GetDiffusionCoefficients () const
 
size_t GetResolution () const
 
real_t GetBoxVolume () const
 
template<typename F >
void AddInitializer (F function)
 
bool IsFixedSubstance ()
 Return if a substance is stationary, e.g. (mu == 0 && dc == 0) More...
 
void SetBoundaryCondition (std::unique_ptr< BoundaryCondition > bc)
 Set the boundary condition, takes ownership of the object. More...
 
BoundaryConditionGetBoundaryCondition () const
 Returns the boundary condition. Does not transfer ownership. More...
 
void SetBoundaryConditionType (BoundaryConditionType bc_type)
 Sets boundary condition type. More...
 
BoundaryConditionType GetBoundaryConditionType () const
 Returns the BoundaryConditionType, see BoundaryConditionType More...
 
void PrintInfo (std::ostream &out=std::cout)
 Print information about the Diffusion Grid. More...
 
void PrintInfoWithInitialization ()
 Print the information after initialization. More...
 
bool IsInitialized () const
 Returns if the grid has been initialized. More...
 
void TurnOffGradientCalculation ()
 
- Public Member Functions inherited from bdm::ScalarField
 ScalarField ()=default
 
 ScalarField (const TRootIOCtor *)
 
 ~ScalarField () override=default
 
 BDM_CLASS_DEF_OVERRIDE (ScalarField, 1)
 
- Public Member Functions inherited from bdm::Continuum
 Continuum ()=default
 
 Continuum (const TRootIOCtor *)
 
virtual ~Continuum ()=default
 
void IntegrateTimeAsynchronously (real_t dt)
 
int GetContinuumId () const
 Returns the ID of the continuum. More...
 
void SetContinuumId (int id)
 Sets the ID of the continuum. More...
 
const std::string & GetContinuumName () const
 Returns the name of the continuum. More...
 
void SetContinuumName (const std::string &name)
 Sets the name of the continuum. More...
 
real_t GetSimulatedTime () const
 Returns the time simulated by the continuum. More...
 
void SetTimeStep (real_t dt)
 
real_t GetTimeStep () const
 Returns the time step for the continuum. More...
 

Private Member Functions

void ParametersCheck (real_t dt)
 
void CopyOldData (const ParallelResizeVector< real_t > &old_c1, const ParallelResizeVector< Real3 > &old_gradients, size_t old_resolution)
 
 BDM_CLASS_DEF_OVERRIDE (DiffusionGrid, 1)
 

Private Attributes

real_t box_length_ = 0
 The side length of each box. More...
 
real_t box_volume_ = 0
 the volume of each box More...
 
ParallelResizeVector< Spinlocklocks_ = {}
 
ParallelResizeVector< real_tc1_ = {}
 The array of concentration values. More...
 
ParallelResizeVector< real_tc2_ = {}
 An extra concentration data buffer for faster value updating. More...
 
ParallelResizeVector< Real3gradients_ = {}
 The array of gradients (x, y, z) More...
 
real_t upper_threshold_ = 1e15
 The maximum concentration value that a box can have. More...
 
real_t lower_threshold_ = 0.0
 The minimum concentration value that a box can have. More...
 
std::array< real_t, 7 > dc_ = {{0}}
 The diffusion coefficients [cc, cw, ce, cs, cn, cb, ct]. More...
 
real_t mu_ = 0
 The decay constant. More...
 
std::array< int32_t, 2 > grid_dimensions_ = {{0}}
 The grid dimensions of the diffusion grid (cubic shaped) More...
 
uint64_t num_boxes_axis_ = 0
 The number of boxes at each axis [x, y, z] (same along each axis) More...
 
size_t total_num_boxes_ = 0
 The total number of boxes in the diffusion grid. More...
 
size_t resolution_ = 0
 
real_t last_dt_ = 0.0
 The last timestep dt used for the diffusion grid update Diffuse(dt) More...
 
bool parity_ = false
 If false, grid dimensions are even; if true, they are odd. More...
 
std::vector< std::function< real_t(real_t, real_t, real_t)> > initializers_
 
bool init_gradient_ = false
 
BoundaryConditionType bc_type_ = BoundaryConditionType::kDirichlet
 Type of boundary conditions. More...
 
std::unique_ptr< BoundaryConditionboundary_condition_ = nullptr
 Object that implements the boundary conditions. More...
 
bool initialized_ = false
 Flag to indicate if the grid is initialized. More...
 
bool print_info_with_initialization_ = false
 
bool precompute_gradients_ = true
 

Friends

class EulerGrid
 
class EulerDepletionGrid
 
class TestGrid
 

Detailed Description

Definition at line 90 of file diffusion_grid.h.

Constructor & Destructor Documentation

◆ DiffusionGrid() [1/5]

bdm::DiffusionGrid::DiffusionGrid ( )
default

◆ DiffusionGrid() [2/5]

bdm::DiffusionGrid::DiffusionGrid ( const TRootIOCtor *  )
inlineexplicit

Definition at line 93 of file diffusion_grid.h.

◆ DiffusionGrid() [3/5]

bdm::DiffusionGrid::DiffusionGrid ( int  substance_id,
const std::string &  substance_name,
real_t  dc,
real_t  mu,
int  resolution = 10 
)

Definition at line 61 of file diffusion_grid.cc.

◆ ~DiffusionGrid()

bdm::DiffusionGrid::~DiffusionGrid ( )
overridedefault

◆ DiffusionGrid() [4/5]

bdm::DiffusionGrid::DiffusionGrid ( const DiffusionGrid )
delete

◆ DiffusionGrid() [5/5]

bdm::DiffusionGrid::DiffusionGrid ( DiffusionGrid &&  )
delete

Member Function Documentation

◆ AddInitializer()

template<typename F >
void bdm::DiffusionGrid::AddInitializer ( function)
inline

Definition at line 311 of file diffusion_grid.h.

◆ BDM_CLASS_DEF_OVERRIDE()

bdm::DiffusionGrid::BDM_CLASS_DEF_OVERRIDE ( DiffusionGrid  ,
 
)
private

◆ CalculateGradient()

void bdm::DiffusionGrid::CalculateGradient ( )

Calculates the gradient for each box in the diffusion grid. The gradient is calculated in each direction (x, y, z) as following:

c(x + box_length_) - c(x - box_length) / (2 * box_length_),

where c(x) implies the concentration at position x

At the edges the gradient is the same as the box next to it

Definition at line 317 of file diffusion_grid.cc.

◆ ChangeConcentrationBy() [1/2]

void bdm::DiffusionGrid::ChangeConcentrationBy ( const Real3 position,
real_t  amount,
InteractionMode  mode = InteractionMode::kAdditive,
bool  scale_with_resolution = false 
)

Increase the concentration at a specified position a with specified amount. The InteractionMode defines how the new concentration is calculated. If the interaction mode is kAdditive, the new concentration is the sum of the old concentration and the amount. If the interaction mode is kExponential, the new concentration is the product of the old concentration and the amount. If the interaction mode is kLogistic, the new concentration is the logistic function of the old concentration and the amount, e.g. additive but scaled with (upper_threshold_ - u) or (u - lower_threshold_) if amount is positive or negative, respectively. Should be used with thresholds 1.0 and 0.0. The parameter scale_with_resolution scales the amount with the inverse of the voxel volume. This is useful if the amount is given in units of, for instance, [mg] but the diffusion grid is in units of [mg / um^3]. It's helpful to model this way to obtain similar/identical results independent of the resolution of the diffusion grid.

Definition at line 355 of file diffusion_grid.cc.

◆ ChangeConcentrationBy() [2/2]

void bdm::DiffusionGrid::ChangeConcentrationBy ( size_t  idx,
real_t  amount,
InteractionMode  mode = InteractionMode::kAdditive,
bool  scale_with_resolution = false 
)

See ChangeConcentrationBy(const Real3&, real_t, InteractionMode, bool).

Increase the concentration at specified box with specified amount.

Definition at line 363 of file diffusion_grid.cc.

◆ CopyOldData()

void bdm::DiffusionGrid::CopyOldData ( const ParallelResizeVector< real_t > &  old_c1,
const ParallelResizeVector< Real3 > &  old_gradients,
size_t  old_resolution 
)
private

Copies the concentration and gradients values to the new (larger) grid. In the 2D case it looks like the following:

                        [0 0  0  0]
          [v1 v2]  -->  [0 v1 v2 0]
          [v3 v4]  -->  [0 v3 v4 0]
                        [0 0  0  0]

The dimensions are doubled in this case from 2x2 to 4x4 If the dimensions would be increased from 2x2 to 3x3, it will still be increased to 4x4 in order for GetBoxIndex to function correctly

Definition at line 211 of file diffusion_grid.cc.

◆ Diffuse()

void bdm::DiffusionGrid::Diffuse ( real_t  dt)

Definition at line 125 of file diffusion_grid.cc.

◆ DiffuseWithClosedEdge()

virtual void bdm::DiffusionGrid::DiffuseWithClosedEdge ( real_t  dt)
pure virtual

◆ DiffuseWithDirichlet()

virtual void bdm::DiffusionGrid::DiffuseWithDirichlet ( real_t  dt)
pure virtual

Compute the diffusion of the substance with Dirichlet boundary conditions ( u = const on the boundary) for a given time step dt.

Implemented in bdm::EulerDepletionGrid, and bdm::EulerGrid.

◆ DiffuseWithNeumann()

virtual void bdm::DiffusionGrid::DiffuseWithNeumann ( real_t  dt)
pure virtual

Compute the diffusion of the substance with Neumann boundary conditions (du/dn = const on the boundary) for a given time step dt.

Implemented in bdm::EulerDepletionGrid, and bdm::EulerGrid.

◆ DiffuseWithOpenEdge()

virtual void bdm::DiffusionGrid::DiffuseWithOpenEdge ( real_t  dt)
pure virtual

◆ DiffuseWithPeriodic()

virtual void bdm::DiffusionGrid::DiffuseWithPeriodic ( real_t  dt)
pure virtual

Compute the diffusion of the substance with Periodic boundary conditions (u_{-1} = {u_N-1}, u_{N} = u_{0}) for a given time step dt.

Implemented in bdm::EulerDepletionGrid, and bdm::EulerGrid.

◆ GetAllConcentrations()

const real_t* bdm::DiffusionGrid::GetAllConcentrations ( ) const
inline

Definition at line 256 of file diffusion_grid.h.

◆ GetAllGradients()

const real_t* bdm::DiffusionGrid::GetAllGradients ( ) const
inline

Definition at line 258 of file diffusion_grid.h.

◆ GetBoundaryCondition()

BoundaryCondition * bdm::DiffusionGrid::GetBoundaryCondition ( ) const

Returns the boundary condition. Does not transfer ownership.

Returns
const Pointer to the boundary condition

Definition at line 536 of file diffusion_grid.cc.

◆ GetBoundaryConditionType()

BoundaryConditionType bdm::DiffusionGrid::GetBoundaryConditionType ( ) const

Returns the BoundaryConditionType, see BoundaryConditionType

Definition at line 548 of file diffusion_grid.cc.

◆ GetBoxCoordinates() [1/2]

std::array< uint32_t, 3 > bdm::DiffusionGrid::GetBoxCoordinates ( const Real3 position) const

Get the coordinates of the box at the specified position.

Definition at line 463 of file diffusion_grid.cc.

◆ GetBoxCoordinates() [2/2]

std::array< uint32_t, 3 > bdm::DiffusionGrid::GetBoxCoordinates ( const size_t  idx) const

Definition at line 485 of file diffusion_grid.cc.

◆ GetBoxIndex() [1/2]

size_t bdm::DiffusionGrid::GetBoxIndex ( const Real3 position) const

Calculates the box index of the substance at specified position.

Definition at line 505 of file diffusion_grid.cc.

◆ GetBoxIndex() [2/2]

size_t bdm::DiffusionGrid::GetBoxIndex ( const std::array< uint32_t, 3 > &  box_coord) const

Calculates the box index at specified box coordinates.

Definition at line 497 of file diffusion_grid.cc.

◆ GetBoxLength()

real_t bdm::DiffusionGrid::GetBoxLength ( ) const
inline

Definition at line 270 of file diffusion_grid.h.

◆ GetBoxVolume()

real_t bdm::DiffusionGrid::GetBoxVolume ( ) const
inline

Definition at line 308 of file diffusion_grid.h.

◆ GetConcentration() [1/2]

real_t bdm::DiffusionGrid::GetConcentration ( const Real3 position) const
inline

Definition at line 170 of file diffusion_grid.h.

◆ GetConcentration() [2/2]

real_t bdm::DiffusionGrid::GetConcentration ( const size_t  idx) const

Get the concentration at specified index.

Get the concentration at specified voxel.

Parameters
idxFlat index of the grid
Returns
c1_[idx]

Definition at line 407 of file diffusion_grid.cc.

◆ GetDecayConstant()

real_t bdm::DiffusionGrid::GetDecayConstant ( ) const
inline

Definition at line 281 of file diffusion_grid.h.

◆ GetDiffusionCoefficients()

const std::array<real_t, 7>& bdm::DiffusionGrid::GetDiffusionCoefficients ( ) const
inline

Definition at line 304 of file diffusion_grid.h.

◆ GetDimensions()

std::array<int32_t, 6> bdm::DiffusionGrid::GetDimensions ( ) const
inline

Definition at line 285 of file diffusion_grid.h.

◆ GetDimensionsPtr()

const int32_t* bdm::DiffusionGrid::GetDimensionsPtr ( ) const
inline

Definition at line 283 of file diffusion_grid.h.

◆ GetGradient() [1/2]

Real3 bdm::DiffusionGrid::GetGradient ( const Real3 position) const
inlineoverridevirtual

Get the gradient at a specified position. By default, the obtained gradient is scaled to norm 1, but with normalize = false one can obtain the full gradient information (e.g. the un-normalized gradient). If the gradient is zero and normalize = true, this method returns a zero vector. The gradient is computed via a central difference scheme on the underlying spatial discretization (central difference for inner boxes, forward/backward difference for edge boxes). The gradients can either be computed once for all boxes (set Param::calculate_gradients to true) or on the fly (set Param::calculate_gradients to false). Depending on the use case, one of the two options might be more efficient (e.g. if the gradient is only needed for a few boxes, 'on the fly' might be more efficient). Please note that the on the fly evaluation of the gradient is done on the current state of the diffusion grid (c1_), thus, it interferes with ChangeConcentrationBy() (e.g. Adding to the substance will change the gradient: evaluating the gradient during the same iteration may yield different results). Further, if you need to visualize the gradient with Paraview, you must set Param::calculate_gradients to true.

Implements bdm::ScalarField.

Definition at line 196 of file diffusion_grid.h.

◆ GetGradient() [2/2]

void bdm::DiffusionGrid::GetGradient ( const Real3 position,
Real3 gradient,
bool  normalize = true 
) const
virtual

Get the gradient at a specified position. See GetGradient(Real3) for details. Note that the pointer to the gradient is passed as an argument and is dereferenced in the method. Thus the gradient must be initialized before calling this method and must not be nullptr.

Definition at line 419 of file diffusion_grid.cc.

◆ GetGridSize()

std::array<int32_t, 3> bdm::DiffusionGrid::GetGridSize ( ) const
inline

Definition at line 296 of file diffusion_grid.h.

◆ GetLastTimestep()

real_t bdm::DiffusionGrid::GetLastTimestep ( ) const
inline

Return the last timestep dt that was used to run Diffuse(dt)

Definition at line 242 of file diffusion_grid.h.

◆ GetLowerThreshold()

real_t bdm::DiffusionGrid::GetLowerThreshold ( ) const
inline

Definition at line 254 of file diffusion_grid.h.

◆ GetNeighboringBoxes() [1/2]

std::array< size_t, 6 > bdm::DiffusionGrid::GetNeighboringBoxes ( size_t  index) const

Determines the indices of the neighboring boxes of a given box index. Order: x-, x+, y-, y+, z-, z+. If the box is on the edge of the grid, the index itself is returned. E.g. if the box is on the left edge of the grid, x- is the same as index.

Definition at line 510 of file diffusion_grid.cc.

◆ GetNeighboringBoxes() [2/2]

std::array< size_t, 6 > bdm::DiffusionGrid::GetNeighboringBoxes ( size_t  index,
const std::array< uint32_t, 3 > &  box_coord 
) const

Determines the indices of the neighboring boxes of given box coordinates. Order: x-, x+, y-, y+, z-, z+. If the box is on the edge of the grid, the index itself is returned. E.g. if the box is on the left edge of the grid, x- is the same as index.

Definition at line 515 of file diffusion_grid.cc.

◆ GetNumBoxes()

size_t bdm::DiffusionGrid::GetNumBoxes ( ) const
inline

Definition at line 268 of file diffusion_grid.h.

◆ GetNumBoxesArray()

std::array<size_t, 3> bdm::DiffusionGrid::GetNumBoxesArray ( ) const
inline

Definition at line 260 of file diffusion_grid.h.

◆ GetResolution()

size_t bdm::DiffusionGrid::GetResolution ( ) const
inline

Definition at line 306 of file diffusion_grid.h.

◆ GetSubstanceId()

int bdm::DiffusionGrid::GetSubstanceId ( ) const
inline

Definition at line 272 of file diffusion_grid.h.

◆ GetSubstanceName()

const std::string& bdm::DiffusionGrid::GetSubstanceName ( ) const
inline

Definition at line 277 of file diffusion_grid.h.

◆ GetUpperThreshold()

real_t bdm::DiffusionGrid::GetUpperThreshold ( ) const
inline

Definition at line 248 of file diffusion_grid.h.

◆ GetValue()

real_t bdm::DiffusionGrid::GetValue ( const Real3 position) const
overridevirtual

Get the value of the scalar field at specified position.

Get the concentration at specified position.

Parameters
position3D position of
Returns
c1_[idx[position]]

Implements bdm::ScalarField.

Definition at line 401 of file diffusion_grid.cc.

◆ Initialize()

void bdm::DiffusionGrid::Initialize ( )
overridevirtual

Initializes the continuum. This method is called via Scheduler::Initialize. For some implementations, this method may be useful, other may not require it. A possibly use case is that agents move in a stationary continuum that is the solution to a timeindependent PDE. In this case, the Initialize method could be used to solve the PDE once at the beginning of the simulation.

Implements bdm::Continuum.

Definition at line 77 of file diffusion_grid.cc.

◆ IsFixedSubstance()

bool bdm::DiffusionGrid::IsFixedSubstance ( )
inline

Return if a substance is stationary, e.g. (mu == 0 && dc == 0)

Definition at line 316 of file diffusion_grid.h.

◆ IsInitialized()

bool bdm::DiffusionGrid::IsInitialized ( ) const
inline

Returns if the grid has been initialized.

Definition at line 345 of file diffusion_grid.h.

◆ operator=() [1/2]

DiffusionGrid& bdm::DiffusionGrid::operator= ( const DiffusionGrid )
delete

◆ operator=() [2/2]

DiffusionGrid& bdm::DiffusionGrid::operator= ( DiffusionGrid &&  )
delete

◆ ParametersCheck()

void bdm::DiffusionGrid::ParametersCheck ( real_t  dt)
private

Definition at line 593 of file diffusion_grid.cc.

◆ PrintInfo()

void bdm::DiffusionGrid::PrintInfo ( std::ostream &  out = std::cout)

Print information about the Diffusion Grid.

Definition at line 552 of file diffusion_grid.cc.

◆ PrintInfoWithInitialization()

void bdm::DiffusionGrid::PrintInfoWithInitialization ( )
inline

Print the information after initialization.

Definition at line 340 of file diffusion_grid.h.

◆ RunInitializers()

void bdm::DiffusionGrid::RunInitializers ( )

Initialize the diffusion grid according to the initialization functions. Note that if your initializers our outside the defined bounds (lower and upper threshold), they will be clamped to the bounds.

Definition at line 253 of file diffusion_grid.cc.

◆ SetBoundaryCondition()

void bdm::DiffusionGrid::SetBoundaryCondition ( std::unique_ptr< BoundaryCondition bc)

Set the boundary condition, takes ownership of the object.

Parameters
bcobject that implements the boundary condition, see for example ConstantBoundaryCondition

Definition at line 531 of file diffusion_grid.cc.

◆ SetBoundaryConditionType()

void bdm::DiffusionGrid::SetBoundaryConditionType ( BoundaryConditionType  bc_type)

Sets boundary condition type.

Definition at line 540 of file diffusion_grid.cc.

◆ SetDecayConstant()

void bdm::DiffusionGrid::SetDecayConstant ( real_t  mu)
inline

Definition at line 236 of file diffusion_grid.h.

◆ SetLowerThreshold()

void bdm::DiffusionGrid::SetLowerThreshold ( real_t  t)
inline

Definition at line 251 of file diffusion_grid.h.

◆ SetUpperThreshold()

void bdm::DiffusionGrid::SetUpperThreshold ( real_t  t)
inline

Definition at line 245 of file diffusion_grid.h.

◆ Step()

void bdm::DiffusionGrid::Step ( real_t  dt)
inlineoverridevirtual

This method integrates the continuum model in time by dt. Typically, the continuum is the solution of a partial differential equation (PDE). Thus, the step method is used to advance the solution in time. If your numerical scheme has requirements for it's parameters depending on dt, make sure to verify them in this method.

Implements bdm::Continuum.

Definition at line 110 of file diffusion_grid.h.

◆ TurnOffGradientCalculation()

void bdm::DiffusionGrid::TurnOffGradientCalculation ( )
inline

Turn off the gradient calculation. Gradients are not precomputed but can be calculated on the fly.

Definition at line 349 of file diffusion_grid.h.

◆ Update()

void bdm::DiffusionGrid::Update ( )
overridevirtual

Updates the grid dimensions, based on the given threshold values. The diffusion grid dimensions need always be larger than the neighbor grid dimensions, so that each simulation object can obtain its local concentration / gradient

Implements bdm::Continuum.

Definition at line 156 of file diffusion_grid.cc.

Friends And Related Function Documentation

◆ EulerDepletionGrid

friend class EulerDepletionGrid
friend

Definition at line 353 of file diffusion_grid.h.

◆ EulerGrid

friend class EulerGrid
friend

Definition at line 352 of file diffusion_grid.h.

◆ TestGrid

friend class TestGrid
friend

Definition at line 354 of file diffusion_grid.h.

Member Data Documentation

◆ bc_type_

BoundaryConditionType bdm::DiffusionGrid::bc_type_ = BoundaryConditionType::kDirichlet
private

Type of boundary conditions.

Definition at line 415 of file diffusion_grid.h.

◆ boundary_condition_

std::unique_ptr<BoundaryCondition> bdm::DiffusionGrid::boundary_condition_ = nullptr
private

Object that implements the boundary conditions.

Definition at line 417 of file diffusion_grid.h.

◆ box_length_

real_t bdm::DiffusionGrid::box_length_ = 0
private

The side length of each box.

Definition at line 375 of file diffusion_grid.h.

◆ box_volume_

real_t bdm::DiffusionGrid::box_volume_ = 0
private

the volume of each box

Definition at line 377 of file diffusion_grid.h.

◆ c1_

ParallelResizeVector<real_t> bdm::DiffusionGrid::c1_ = {}
private

The array of concentration values.

Definition at line 382 of file diffusion_grid.h.

◆ c2_

ParallelResizeVector<real_t> bdm::DiffusionGrid::c2_ = {}
private

An extra concentration data buffer for faster value updating.

Definition at line 384 of file diffusion_grid.h.

◆ dc_

std::array<real_t, 7> bdm::DiffusionGrid::dc_ = {{0}}
private

The diffusion coefficients [cc, cw, ce, cs, cn, cb, ct].

Definition at line 392 of file diffusion_grid.h.

◆ gradients_

ParallelResizeVector<Real3> bdm::DiffusionGrid::gradients_ = {}
private

The array of gradients (x, y, z)

Definition at line 386 of file diffusion_grid.h.

◆ grid_dimensions_

std::array<int32_t, 2> bdm::DiffusionGrid::grid_dimensions_ = {{0}}
private

The grid dimensions of the diffusion grid (cubic shaped)

Definition at line 396 of file diffusion_grid.h.

◆ init_gradient_

bool bdm::DiffusionGrid::init_gradient_ = false
private

Definition at line 413 of file diffusion_grid.h.

◆ initialized_

bool bdm::DiffusionGrid::initialized_ = false
private

Flag to indicate if the grid is initialized.

Definition at line 419 of file diffusion_grid.h.

◆ initializers_

std::vector<std::function<real_t(real_t, real_t, real_t)> > bdm::DiffusionGrid::initializers_
private
Initial value:
=
{}

A list of functions that initialize this diffusion grid ROOT currently doesn't support IO of std::function

Definition at line 410 of file diffusion_grid.h.

◆ last_dt_

real_t bdm::DiffusionGrid::last_dt_ = 0.0
private

The last timestep dt used for the diffusion grid update Diffuse(dt)

Definition at line 405 of file diffusion_grid.h.

◆ locks_

ParallelResizeVector<Spinlock> bdm::DiffusionGrid::locks_ = {}
mutableprivate

Lock for each voxel used to prevent race conditions between multiple threads

Definition at line 380 of file diffusion_grid.h.

◆ lower_threshold_

real_t bdm::DiffusionGrid::lower_threshold_ = 0.0
private

The minimum concentration value that a box can have.

Definition at line 390 of file diffusion_grid.h.

◆ mu_

real_t bdm::DiffusionGrid::mu_ = 0
private

The decay constant.

Definition at line 394 of file diffusion_grid.h.

◆ num_boxes_axis_

uint64_t bdm::DiffusionGrid::num_boxes_axis_ = 0
private

The number of boxes at each axis [x, y, z] (same along each axis)

Definition at line 398 of file diffusion_grid.h.

◆ parity_

bool bdm::DiffusionGrid::parity_ = false
private

If false, grid dimensions are even; if true, they are odd.

Definition at line 407 of file diffusion_grid.h.

◆ precompute_gradients_

bool bdm::DiffusionGrid::precompute_gradients_ = true
private

Flag to avoid gradient computation if not needed. (E.g. if multiple DGs are used but the gradient is only needed for one of them.)

Definition at line 425 of file diffusion_grid.h.

◆ print_info_with_initialization_

bool bdm::DiffusionGrid::print_info_with_initialization_ = false
private

Flag to indicate if we want to print information about the grid after initialization

Definition at line 422 of file diffusion_grid.h.

◆ resolution_

size_t bdm::DiffusionGrid::resolution_ = 0
private

The resolution of the diffusion grid (i.e. number of boxes along each axis)

Definition at line 403 of file diffusion_grid.h.

◆ total_num_boxes_

size_t bdm::DiffusionGrid::total_num_boxes_ = 0
private

The total number of boxes in the diffusion grid.

Definition at line 400 of file diffusion_grid.h.

◆ upper_threshold_

real_t bdm::DiffusionGrid::upper_threshold_ = 1e15
private

The maximum concentration value that a box can have.

Definition at line 388 of file diffusion_grid.h.


The documentation for this class was generated from the following files: