solidificationTemplates.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) 2017 OpenFOAM Foundation
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 "geometricOneField.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class AlphaFieldType, class RhoFieldType>
34 void Foam::porosityModels::solidification::apply
35 (
36  scalarField& Udiag,
37  const scalarField& V,
38  const AlphaFieldType& alpha,
39  const RhoFieldType& rho,
40  const volVectorField& U
41 ) const
42 {
43  const auto& T = mesh_.lookupObject<volScalarField>
44  (
45  IOobject::groupName(TName_, U.group())
46  );
47 
48  for (const label zonei : cellZoneIDs_)
49  {
50  const labelList& cells = mesh_.cellZones()[zonei];
51 
52  for (const label celli : cells)
53  {
54  Udiag[celli] +=
55  V[celli]*alpha[celli]*rho[celli]*D_->value(T[celli]);
56  }
57  }
58 }
59 
60 
61 template<class AlphaFieldType, class RhoFieldType>
62 void Foam::porosityModels::solidification::apply
63 (
64  tensorField& AU,
65  const AlphaFieldType& alpha,
66  const RhoFieldType& rho,
67  const volVectorField& U
68 ) const
69 {
70  const auto& T = mesh_.lookupObject<volScalarField>
71  (
72  IOobject::groupName(TName_, U.group())
73  );
74 
75  for (const label zonei : cellZoneIDs_)
76  {
77  const labelList& cells = mesh_.cellZones()[zonei];
78 
79  for (const label celli : cells)
80  {
81  AU[celli] +=
82  tensor::I*alpha[celli]*rho[celli]*D_->value(T[celli]);
83  }
84  }
85 }
86 
87 
88 template<class RhoFieldType>
89 void Foam::porosityModels::solidification::apply
90 (
91  scalarField& Udiag,
92  const scalarField& V,
93  const RhoFieldType& rho,
94  const volVectorField& U
95 ) const
96 {
97  if (alphaName_ == "none")
98  {
99  return apply(Udiag, V, geometricOneField(), rho, U);
100  }
101  else
102  {
103  const auto& alpha = mesh_.lookupObject<volScalarField>
104  (
105  IOobject::groupName(alphaName_, U.group())
106  );
107 
108  return apply(Udiag, V, alpha, rho, U);
109  }
110 }
111 
112 
113 template<class RhoFieldType>
114 void Foam::porosityModels::solidification::apply
115 (
116  tensorField& AU,
117  const RhoFieldType& rho,
118  const volVectorField& U
119 ) const
120 {
121  if (alphaName_ == "none")
122  {
123  return apply(AU, geometricOneField(), rho, U);
124  }
125  else
126  {
127  const auto& alpha = mesh_.lookupObject<volScalarField>
128  (
129  IOobject::groupName(alphaName_, U.group())
130  );
131 
132  return apply(AU, alpha, rho, U);
133  }
134 }
135 
136 
137 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Tensor::I
static const Tensor I
Definition: Tensor.H:85
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::porosityModel::cellZoneIDs_
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:93
rho
rho
Definition: readInitialConditions.H:88
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
geometricOneField.H
Foam::apply
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Definition: parcelSelectionDetail.C:82
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::porosityModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:78