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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "fvMesh.H"
30#include "surfaceFields.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
40 (
44 );
45
47 (
50 geometry
51 );
52}
53
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
58(
59 const lduMesh& mesh,
61)
62:
64{
65 const fvMesh& fvmesh = refCast<const fvMesh>(mesh);
66
67 //agglomerate(mesh, sqrt(fvmesh.magSf().primitiveField()));
69 (
70 mesh,
71 mag
72 (
74 (
75 fvmesh.Sf().primitiveField()
76 /sqrt(fvmesh.magSf().primitiveField()),
77 vector(1, 1.01, 1.02)
78 //vector::one
79 )
80 )
81 );
82}
83
84
86(
87 const lduMesh& mesh,
88 const scalarField& cellVolumes,
89 const vectorField& faceAreas,
91)
92:
94{
95 //agglomerate(mesh, sqrt(mag(faceAreas)));
97 (
98 mesh,
99 mag
100 (
102 (
103 faceAreas
104 /sqrt(mag(faceAreas)),
105 vector(1, 1.01, 1.02)
106 //vector::one
107 )
108 )
109 );
110}
111
112
113// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Geometric agglomerated algebraic multigrid agglomeration class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Agglomerate using the pair algorithm.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
Agglomerate using the pair algorithm.
void agglomerate(const lduMesh &mesh, const scalarField &faceWeights)
Agglomerate all levels starting from the given face weights.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
runTime controlDict().readEntry("adjustTimeStep"
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Vector< scalar > vector
Definition: vector.H:61
Foam::surfaceFields.