faceAreaPairGAMGAgglomeration.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvMesh.H"
31 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(faceAreaPairGAMGAgglomeration, 0);
39 
41  (
42  GAMGAgglomeration,
43  faceAreaPairGAMGAgglomeration,
44  lduMesh
45  );
46 
48  (
49  GAMGAgglomeration,
50  faceAreaPairGAMGAgglomeration,
51  geometry
52  );
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const lduMesh& mesh,
61  const dictionary& controlDict
62 )
63 :
65 {
66  const fvMesh& fvmesh = refCast<const fvMesh>(mesh);
67 
68  //agglomerate(mesh, sqrt(fvmesh.magSf().primitiveField()));
70  (
72  0, //mesh,
73  mag
74  (
76  (
77  fvmesh.Sf().primitiveField()
78  /sqrt(fvmesh.magSf().primitiveField()),
79  vector(1, 1.01, 1.02)
80  //vector::one
81  )
82  ),
83  true
84  );
85 }
86 
87 
89 (
90  const lduMesh& mesh,
91  const scalarField& cellVolumes,
92  const vectorField& faceAreas,
93  const dictionary& controlDict
94 )
95 :
97 {
98  //agglomerate(mesh, sqrt(mag(faceAreas)));
100  (
102  0, //mesh,
103  mag
104  (
106  (
107  faceAreas
108  /sqrt(mag(faceAreas)),
109  vector(1, 1.01, 1.02)
110  //vector::one
111  )
112  ),
113  true
114  );
115 }
116 
117 
118 // ************************************************************************* //
Foam::surfaceFields.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Agglomerate using the pair algorithm.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
Macros for easy insertion into run-time selection tables.
faceAreaPairGAMGAgglomeration(const lduMesh &mesh, const dictionary &controlDict)
Construct given mesh and controls.
static tmp< labelField > agglomerate(label &nCoarseCells, const lduAddressing &fineMatrixAddressing, const scalarField &faceWeights)
Calculate and return agglomeration.
dynamicFvMesh & mesh
const lduMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:157
Vector< scalar > vector
Definition: vector.H:57
runTime controlDict().readEntry("adjustTimeStep"
The central control dictionary, the contents of which are either taken directly from the FOAM_CONTROL...
Definition: debug.C:142
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
defineTypeNameAndDebug(combustionModel, 0)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label nCellsInCoarsestLevel_
Number of cells in coarsest level.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)