OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

porousZone Class Reference

Collaboration diagram for porousZone:

Collaboration graph
[legend]
List of all members.

Detailed Description

Porous zone definition based on cell zones.

Porous zone definition based on cell zones and parameters obtained from a control dictionary constructed from the given stream. The orientation of the porous region is defined with the same notation as a coordinateSystem, but only a Cartesian coordinate system is valid.

Implemented porosity models:

powerLaw (C0 and C1 parameters)

\[ S = - \rho C_0 |U|^{(C_1 - 1)/2} U \]

Darcy-Forchheimer (d and f parameters)

\[ S = - (\mu \, d \, U + \frac{\rho |U|}{2} \, f) U \]

Since negative Darcy/Forchheimer parameters are invalid, they can be used to specify a multiplier (of the max component).

The porousZones method porousZones::ddt() mirrors the normal fvm::ddt() method, but accounts for the effective volume of the cells.

See also:
porousZones and coordinateSystems
Source files

Definition at line 85 of file porousZone.H.


Public Member Functions

 porousZone (const fvMesh &, const word &name, const dictionary &)
 Construct from components.
autoPtr< porousZoneclone () const
 Return clone.
const wordzoneName () const
 cellZone name
label zoneId () const
 cellZone number
const dictionarydict () const
 dictionary values used for the porousZone
const coordinateSystemcoordSys () const
 Return coordinate system.
const pointorigin () const
 Return origin.
const vectoraxis () const
 Return axis.
const scalarporosity () const
 Return porosity.
scalarporosity ()
 Edit access to porosity.
template<class Type>
void modifyDdt (fvMatrix< Type > &) const
 modify time derivative elements according to porosity
void addResistance (fvVectorMatrix &UEqn) const
 Add the viscous and inertial resistance force contribution.
void addResistance (const fvVectorMatrix &UEqn, volTensorField &AU, bool correctAUprocBC=true) const
 Add the viscous and inertial resistance force contribution.
void writeDict (Ostream &, bool subDict=true) const
 Write the porousZone dictionary.

Friends

Ostreamoperator<< (Ostream &, const porousZone &)

Classes

class  iNew
 Return a pointer to a new porousZone created on freestore. More...

Constructor & Destructor Documentation

porousZone ( const fvMesh ,
const word name,
const dictionary  
)

Construct from components.

Definition at line 56 of file porousZone.C.


Member Function Documentation

autoPtr<porousZone> clone (  )  const [inline]

Return clone.

Definition at line 191 of file porousZone.H.

const word& zoneName (  )  const [inline]

cellZone name

Definition at line 240 of file porousZone.H.

label zoneId (  )  const [inline]

cellZone number

Definition at line 246 of file porousZone.H.

const dictionary& dict (  )  const [inline]

dictionary values used for the porousZone

Definition at line 252 of file porousZone.H.

const coordinateSystem& coordSys (  )  const [inline]

Return coordinate system.

Definition at line 258 of file porousZone.H.

const point& origin (  )  const [inline]

Return origin.

Definition at line 264 of file porousZone.H.

References coordinateSystem::origin().

Here is the call graph for this function:

const vector& axis (  )  const [inline]

Return axis.

Definition at line 270 of file porousZone.H.

References coordinateSystem::axis().

Here is the call graph for this function:

const scalar& porosity (  )  const [inline]

Return porosity.

Definition at line 276 of file porousZone.H.

scalar& porosity (  )  [inline]

Edit access to porosity.

Definition at line 282 of file porousZone.H.

void modifyDdt ( fvMatrix< Type > &   )  const

modify time derivative elements according to porosity

Definition at line 24 of file porousZoneTemplates.C.

References polyMesh::cellZones(), lduMatrix::diag(), forAll, and fvMatrix::source().

Here is the call graph for this function:

void addResistance ( fvVectorMatrix UEqn  )  const

Add the viscous and inertial resistance force contribution.

to the momentum equation

Definition at line 185 of file porousZone.C.

void addResistance ( const fvVectorMatrix UEqn,
volTensorField AU,
bool  correctAUprocBC = true 
) const

Add the viscous and inertial resistance force contribution.

to the tensorial diagonal. Optionally correct the processor BCs of AU.

Definition at line 266 of file porousZone.C.

void writeDict ( Ostream ,
bool  subDict = true 
) const

Write the porousZone dictionary.

Definition at line 349 of file porousZone.C.


Friends And Related Function Documentation

Ostream& operator<< ( Ostream os,
const porousZone pZone 
) [friend]

Definition at line 393 of file porousZone.C.


The documentation for this class was generated from the following files: Copyright © 2000-2008 OpenCFD Ltd