jouleHeatingSourceTemplates.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) 2016-2024 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 "emptyFvPatchField.H"
29 
30 template<class Type>
31 void Foam::fv::jouleHeatingSource::initialiseSigma
32 (
33  const dictionary& dict,
34  autoPtr<Function1<Type>>& sigmaVsTPtr
35 )
36 {
37  typedef GeometricField<Type, fvPatchField, volMesh> FieldType;
38 
39  IOobject io
40  (
41  IOobject::scopedName(typeName, "sigma"),
42  mesh_.time().timeName(),
43  mesh_.thisDb(),
47  );
48 
49  autoPtr<FieldType> tsigma;
50 
51  if (dict.found("sigma"))
52  {
53  // Sigma to be defined using a Function1 type
54  sigmaVsTPtr = Function1<Type>::New("sigma", dict, &mesh_);
55 
56  tsigma.reset
57  (
58  new FieldType
59  (
60  io,
61  mesh_,
62  Foam::zero{}, // value
64  )
65  );
66 
67  Info<< " Conductivity 'sigma' read from dictionary as f(T)"
68  << nl << endl;
69  }
70  else
71  {
72  // Sigma to be defined by user input
73  io.readOpt(IOobject::MUST_READ);
74 
75  tsigma.reset(new FieldType(io, mesh_));
76 
77  Info<< " Conductivity 'sigma' read from file" << nl << endl;
78  }
79 
81 }
82 
83 
84 template<class Type>
86 Foam::fv::jouleHeatingSource::updateSigma
87 (
88  const autoPtr<Function1<Type>>& sigmaVsTPtr
89 ) const
90 {
92 
93  auto& sigma = mesh_.lookupObjectRef<FieldType>
94  (
95  IOobject::scopedName(typeName, "sigma")
96  );
97 
98  if (!sigmaVsTPtr)
99  {
100  // Electrical conductivity field, sigma, was specified by the user
101  return sigma;
102  }
103 
104  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
105 
106  // Internal field
107  forAll(sigma, i)
108  {
109  sigma[i] = sigmaVsTPtr->value(T[i]);
110  }
111 
112 
113  // Boundary field
114  auto& bf = sigma.boundaryFieldRef();
115  forAll(bf, patchi)
116  {
117  fvPatchField<Type>& pf = bf[patchi];
118  if (!isA<emptyFvPatchField<Type>>(pf))
119  {
120  const scalarField& Tbf = T.boundaryField()[patchi];
121  forAll(pf, facei)
122  {
123  pf[facei] = sigmaVsTPtr->value(Tbf[facei]);
124  }
125  }
126  }
127 
128  // Update processor patches
129  sigma.correctBoundaryConditions();
130 
131  return sigma;
132 }
133 
134 
135 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimPower
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
Definition: dimensionSets.H:54
Nothing to be read.
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition: typeInfo.H:87
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Request registration (bool: true)