solidificationTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "volFields.H"
27 #include "geometricOneField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class AlphaFieldType, class RhoFieldType>
32 void Foam::porosityModels::solidification::apply
33 (
34  scalarField& Udiag,
35  const scalarField& V,
36  const AlphaFieldType& alpha,
37  const RhoFieldType& rho,
38  const volVectorField& U
39 ) const
40 {
42 
44  (
45  IOobject::groupName(TName_, U.group())
46  );
47 
48  forAll(cells, i)
49  {
50  const label celli = cells[i];
51  Udiag[celli] +=
52  V[celli]*alpha[celli]*rho[celli]*D_->value(T[celli]);
53  }
54 }
55 
56 
57 template<class AlphaFieldType, class RhoFieldType>
58 void Foam::porosityModels::solidification::apply
59 (
60  tensorField& AU,
61  const AlphaFieldType& alpha,
62  const RhoFieldType& rho,
63  const volVectorField& U
64 ) const
65 {
66  const labelList& cells = mesh_.cellZones()[zoneName_];
67 
68  const volScalarField& T = mesh_.lookupObject<volScalarField>
69  (
70  IOobject::groupName(TName_, U.group())
71  );
72 
73  forAll(cells, i)
74  {
75  const label celli = cells[i];
76  AU[celli] +=
77  tensor::I*alpha[celli]*rho[celli]*D_->value(T[celli]);
78  }
79 }
80 
81 
82 template<class RhoFieldType>
83 void Foam::porosityModels::solidification::apply
84 (
85  scalarField& Udiag,
86  const scalarField& V,
87  const RhoFieldType& rho,
88  const volVectorField& U
89 ) const
90 {
91  if (alphaName_ == "none")
92  {
93  return apply(Udiag, V, geometricOneField(), rho, U);
94  }
95  else
96  {
97  const volScalarField& alpha = mesh_.lookupObject<volScalarField>
98  (
99  IOobject::groupName(alphaName_, U.group())
100  );
101 
102  return apply(Udiag, V, alpha, rho, U);
103  }
104 }
105 
106 
107 template<class RhoFieldType>
108 void Foam::porosityModels::solidification::apply
109 (
110  tensorField& AU,
111  const RhoFieldType& rho,
112  const volVectorField& U
113 ) const
114 {
115  if (alphaName_ == "none")
116  {
117  return apply(AU, geometricOneField(), rho, U);
118  }
119  else
120  {
121  const volScalarField& alpha = mesh_.lookupObject<volScalarField>
122  (
123  IOobject::groupName(alphaName_, U.group())
124  );
125 
126  return apply(AU, alpha, rho, U);
127  }
128 }
129 
130 
131 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static word groupName(Name name, const word &group)
static const Tensor I
Definition: Tensor.H:83
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:75
word zoneName_
Name of cellZone.
Definition: porosityModel.H:84
volScalarField scalarField(fieldObject, mesh)
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)