53 const word& modelType,
58 zoneName_(
name +
":porousZone")
60 const dictionary& coeffs(
dict.optionalSubDict(modelType +
"Coeffs"));
66 scalar searchSpan(coeffs.get<scalar>(
"searchSpan"));
69 const word topSurfaceFileName(coeffs.get<
word>(
"topSurface"));
79 Function1<scalar>::New(
"Sigma",
dict, &
mesh)
88 mesh.time().constant(),
97 boundBox surfBounds(searchSurf.points());
100 surfBounds.min() - searchSpan*zDir, surfBounds.max()
107 label porousCelli = 0;
111 if (searchBounds.contains(
C[celli]))
113 porousCells[porousCelli++] = celli;
117 porousCells.setSize(porousCelli);
122 forAll(porousCells, porousCelli)
124 start[porousCelli] =
C[porousCells[porousCelli]];
125 end[porousCelli] = start[porousCelli] + searchSpan*zDir;
134 List<pointIndexHit> hitInfo;
135 searchSurf.findLine(start, end, hitInfo);
139 forAll(porousCells, celli)
145 porousCells[porousCelli] = porousCells[celli];
148 (hit.point() -
C[porousCells[porousCelli]]) & zDir;
155 porousCells.setSize(porousCelli);
156 zTop.setSize(porousCelli);
159 triSurface groundSurface
164 mesh.boundaryMesh().patchSet(groundPatches),
182 forAll(groundSurfaceProcTris, i)
184 nTris += groundSurfaceProcTris[i].size();
187 List<labelledTri> groundSurfaceTris(nTris);
190 forAll(groundSurfaceProcTris, i)
192 forAll(groundSurfaceProcTris[i], j)
194 groundSurfaceTris[trii] = groundSurfaceProcTris[i][j];
195 groundSurfaceTris[trii][0] += offset;
196 groundSurfaceTris[trii][1] += offset;
197 groundSurfaceTris[trii][2] += offset;
200 offset += groundSurfaceProcPoints[i].size();
204 forAll(groundSurfaceProcPoints, i)
206 nPoints += groundSurfaceProcPoints[i].size();
212 forAll(groundSurfaceProcPoints, i)
214 forAll(groundSurfaceProcPoints[i], j)
216 groundSurfacePoints[pointi++] = groundSurfaceProcPoints[i][j];
220 groundSurface = triSurface(groundSurfaceTris, groundSurfacePoints);
224 triSurfaceSearch groundSearch(groundSurface);
228 start.setSize(porousCelli);
229 end.setSize(porousCelli);
231 forAll(porousCells, porousCelli)
233 start[porousCelli] =
C[porousCells[porousCelli]];
234 end[porousCelli] = start[porousCelli] - searchSpan*zDir;
237 groundSearch.findLine(start, end, hitInfo);
241 forAll(porousCells, porousCelli)
247 zBottom[porousCelli] =
248 (
C[porousCells[porousCelli]] - hit.point()) & zDir;
256 Sigma_ = SigmaFunc->value(zNorm);
266 zoneID = cellZones.size();
267 cellZones.setSize(zoneID + 1);
285 <<
": zone already exists"
291Foam::porosityModels::powerLawLopesdaCosta::powerLawLopesdaCosta
294 const word& modelType,
297 const word& dummyCellZoneName
307 powerLawLopesdaCostaZone::zoneName_
309 Cd_(coeffs_.get<scalar>(
"Cd")),
413 dict_.writeEntry(name_,
os);
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
bool hit() const noexcept
Is there a hit?
const point_type & point() const noexcept
Return point, no checks.
const Field< point_type > & points() const noexcept
Return reference to global points.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
void setSize(const label n)
Same as resize().
void size(const label n)
Older name for setAddressableSize.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
static bool & parRun() noexcept
Test if this a parallel run.
label size() const noexcept
The number of entries in the list.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bounding box defined in terms of min/max extrema points.
const point & max() const noexcept
Maximum describing the bounding box.
const point & min() const noexcept
Minimum describing the bounding box.
bool contains(const point &pt) const
Contains point? (inside or on edge).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Mesh data needed to do the Finite Volume discretisation.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Top level model for porosity models.
const fvMesh & mesh_
Reference to the mesh database.
dictionary coeffs_
Model coefficients dictionary.
porosityModel(const porosityModel &)=delete
No copy construct.
const dictionary dict_
Dictionary used for model construction.
const dictionary & dict() const
Return dictionary used for model construction.
const word & name() const
Return const access to the porosity model name.
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu)
Return the force over the cell zone(s).
const word zoneName_
Automatically generated zone name for this porous zone.
scalarField Sigma_
Porosity surface area per unit volume zone field.
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
powerLawLopesdaCostaZone(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Constructor.
Variant of the power law porosity model with spatially varying drag coefficient.
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
bool writeData(Ostream &os) const
Write.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
virtual ~powerLawLopesdaCosta()
Destructor.
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
IOoject and searching on triSurface.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual tmp< pointField > points() const
Get the points that define the surface.
Helper class to search on triSurface.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
fvMatrix< vector > fvVectorMatrix
errorManip< error > abort(error &err)
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
#define forAll(list, i)
Loop across all elements in list.