normImpl.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) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
34 bool Foam::functionObjects::norm::calcNorm()
35 {
36  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
37  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
38  typedef DimensionedField<Type, polySurfaceGeoMesh> SurfFieldType;
39 
40  if (mesh_.foundObject<VolFieldType>(fieldName_))
41  {
42  return store
43  (
45  calcNormType<VolFieldType>()
46  );
47  }
48  else if (mesh_.foundObject<SurfaceFieldType>(fieldName_))
49  {
50  return store
51  (
53  calcNormType<SurfaceFieldType>()
54  );
55  }
56  else if (mesh_.foundObject<SurfFieldType>(fieldName_))
57  {
58  return store
59  (
61  calcNormType<SurfFieldType>()
62  );
63  }
64 
65  return false;
66 }
67 
68 
69 template<class GeoFieldType>
70 Foam::tmp<GeoFieldType> Foam::functionObjects::norm::calcNormType()
71 {
72  const GeoFieldType& field = mesh_.lookupObject<GeoFieldType>(fieldName_);
73 
74  const dimensionedScalar perturb(field.dimensions(), SMALL);
75 
76  switch (norm_)
77  {
78  case normType::L1:
79  {
80  return field/stabilise(sumMag(field), perturb);
81  }
82 
83  case normType::L2:
84  {
85  return field/stabilise(mag(field), perturb);
86  }
87 
88  case normType::LP:
89  {
90  return
91  field
92  /stabilise
93  (
94  pow(pow(mag(field), p_), scalar(1)/p_),
95  perturb
96  );
97  }
98 
99  case normType::MAX:
100  {
101  return field/stabilise(max(mag(field)), perturb);
102  }
103 
104  case normType::COMPOSITE:
105  {
106  const scalar t = mesh_.time().timeOutputValue();
107 
108  const dimensionedScalar divisor
109  (
110  field.dimensions(),
111  divisorPtr_->value(t)
112  );
113 
114  return field/stabilise(divisor, perturb);
115  }
116 
117  case normType::FIELD:
118  {
119  return field/stabilise(fieldNorm(field), perturb);
120  }
121 
122  default:
123  break;
124  }
126  return nullptr;
127 }
128 
129 
130 template<class Type>
131 Foam::tmp<Foam::volScalarField> Foam::functionObjects::norm::fieldNorm
132 (
134 )
135 {
136  return mesh_.lookupObject<volScalarField>(divisorFieldName_);
137 }
138 
139 
140 template<class Type>
141 Foam::tmp<Foam::surfaceScalarField> Foam::functionObjects::norm::fieldNorm
142 (
144 )
145 {
146  return mesh_.lookupObject<surfaceScalarField>(divisorFieldName_);
147 }
148 
149 
150 template<class Type>
151 Foam::tmp<Foam::polySurfaceScalarField> Foam::functionObjects::norm::fieldNorm
152 (
154 )
155 {
156  return mesh_.lookupObject<polySurfaceScalarField>(divisorFieldName_);
157 }
158 
159 
160 // ************************************************************************* //
Foam::surfaceFields.
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1)
rDeltaTY field()
word resultName_
Name of result field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word fieldName_
Name of field to process.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const fvMesh & mesh_
Reference to the fvMesh.