Loading...
Searching...
No Matches
sigmoidalHeaviside.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) 2007-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "sigmoidalHeaviside.H"
32#include "volFields.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37namespace Foam
38{
39
40// * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * * //
41
44(
45 topOInterpolationFunction,
46 sigmoidalHeaviside,
48);
49
50// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
51
53{
54 scalar averageVol(gAverage(mesh_.V().field()));
55 const Vector<label>& solutionD = mesh_.solutionD();
56 const boundBox& bounds = mesh_.bounds();
57 forAll(solutionD, idir)
58 {
59 if (solutionD[idir] == -1)
60 {
61 averageVol /= bounds.span()[idir];
62 break;
63 }
64 }
65
66 scalar width = pow(averageVol, scalar(1)/scalar(mesh_.nGeometricD()));
67 scalar multMeanRadius =
69 (
70 "meanRadiusMult", {{"scale", 2306}}, 2
71 );
73 << "Computed near-band width :: " << width
74 << " and multiplying with " << multMeanRadius << endl;
76 return multMeanRadius*width;
77}
78
79
80// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
81
82sigmoidalHeaviside::sigmoidalHeaviside
83(
84 const fvMesh& mesh,
85 const dictionary& dict
86)
87:
89 dNB_(dict.getOrDefault<scalar>("d", computeNearBandWidth()))
90 //dNB_(Function1<scalar>::New("d", dict))
91{}
92
93
94// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
95
97(
98 const scalarField& arg,
99 scalarField& res
100) const
101{
102 const scalar timeValue = mesh_.time().timeOutputValue();
103 const scalar t(timeValue == 0 ? 1. : timeValue);
104 const scalar pi = constant::mathematical::pi;
106 << type() << "::interpolate:: t, dNB " << t << ", " << dNB_ << endl;
107
108 res = 0.5*(scalar(1) + arg/dNB_ + sin(pi*arg/dNB_)/pi);
109 res = max(min(scalar(1), res), scalar(0));
110}
111
112
114{
116 scalarField& deriv = tderiv.ref();
117 const scalar timeValue = mesh_.time().timeOutputValue();
118 const scalar t(timeValue == 0 ? 1. : timeValue);
119 const scalar pi = constant::mathematical::pi;
120 scalarField argLimited(max(min(dNB_, arg), -dNB_));
122 << type() << "::derivative:: t, dNB " << t << ", " << dNB_ << endl;
123
124 deriv = 0.5*(scalar(1) + cos(pi*argLimited/dNB_))/dNB_;
125
126 return tderiv;
127}
128
129
130// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131
132} // End namespace Foam
133
134// ************************************************************************* //
constexpr scalar pi(M_PI)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion).
Definition TimeStateI.H:24
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value using any compatibility names if needed.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition polyMesh.C:861
A smooth Heaviside function to project the signed distance field in level set topology optimization.
virtual tmp< scalarField > derivative(const scalarField &arg) const
Return of function with respect to the argument field.
virtual void interpolate(const scalarField &arg, scalarField &res) const
Interpolate argument to result.
scalar computeNearBandWidth() const
Compute the near-band width of the fluid-solid interface as.
scalar dNB_
The near-band distance.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
Namespace for bounding specifications. At the moment, mostly for tables.
constexpr scalar pi(M_PI)
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar sin(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299