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-2020 OpenCFD Ltd.
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
28#include "emptyFvPatchField.H"
29
30template<class Type>
31void Foam::fv::jouleHeatingSource::initialiseSigma
32(
33 const dictionary& dict,
34 autoPtr<Function1<Type>>& sigmaVsTPtr
35)
36{
37 typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
38
39 if (dict.found("sigma"))
40 {
41 // Sigma to be defined using a Function1 type
42 sigmaVsTPtr = Function1<Type>::New("sigma", dict, &mesh_);
43
44 auto tsigma = tmp<VolFieldType>::New
45 (
46 IOobject
47 (
48 typeName + ":sigma",
50 mesh_,
53 ),
54 mesh_,
55 dimensioned<Type>(sqr(dimCurrent)/dimPower/dimLength, Zero)
56 );
57
58 mesh_.objectRegistry::store(tsigma.ptr());
59
60 Info<< " Conductivity 'sigma' read from dictionary as f(T)"
61 << nl << endl;
62 }
63 else
64 {
65 // Sigma to be defined by user input
66 auto tsigma = tmp<VolFieldType>::New
67 (
68 IOobject
69 (
70 typeName + ":sigma",
72 mesh_,
75 ),
76 mesh_
77 );
78
79 mesh_.objectRegistry::store(tsigma.ptr());
80
81 Info<< " Conductivity 'sigma' read from file" << nl << endl;
82 }
83}
84
85
86template<class Type>
88Foam::fv::jouleHeatingSource::updateSigma
89(
90 const autoPtr<Function1<Type>>& sigmaVsTPtr
91) const
92{
94
95 auto& sigma = mesh_.lookupObjectRef<VolFieldType>(typeName + ":sigma");
96
97 if (!sigmaVsTPtr)
98 {
99 // Electrical conductivity field, sigma, was specified by the user
100 return sigma;
101 }
102
103 const auto& T = mesh_.lookupObject<volScalarField>(TName_);
104
105 // Internal field
106 forAll(sigma, i)
107 {
108 sigma[i] = sigmaVsTPtr->value(T[i]);
109 }
110
111
112 // Boundary field
113 typename VolFieldType::Boundary& bf = sigma.boundaryFieldRef();
114 forAll(bf, patchi)
115 {
116 fvPatchField<Type>& pf = bf[patchi];
118 {
119 const scalarField& Tbf = T.boundaryField()[patchi];
120 forAll(pf, facei)
121 {
122 pf[facei] = sigmaVsTPtr->value(Tbf[facei]);
123 }
124 }
125 }
126
127 // Update processor patches
128 sigma.correctBoundaryConditions();
129
130 return sigma;
131}
132
133
134// ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
Generic GeometricField class.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const Type & value() const
Return const reference to value.
This boundary condition provides an 'empty' condition for reduced dimensions cases,...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
const volScalarField & T
const dimensionSet dimPower
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
Definition: dimensionSets.H:56
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333