All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2018 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  (
43  IOobject::groupName(TName_, U.group())
44  );
45 
46  forAll(cellZoneIDs_, zoneI)
47  {
48  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zoneI]];
49 
50  forAll(cells, i)
51  {
52  const label celli = cells[i];
53  Udiag[celli] +=
54  V[celli]*alpha[celli]*rho[celli]*D_->value(T[celli]);
55  }
56  }
57 }
58 
59 
60 template<class AlphaFieldType, class RhoFieldType>
61 void Foam::porosityModels::solidification::apply
62 (
63  tensorField& AU,
64  const AlphaFieldType& alpha,
65  const RhoFieldType& rho,
66  const volVectorField& U
67 ) const
68 {
70  (
71  IOobject::groupName(TName_, U.group())
72  );
73 
74  forAll(cellZoneIDs_, zoneI)
75  {
76  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zoneI]];
77 
78  forAll(cells, i)
79  {
80  const label celli = cells[i];
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  {
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  {
128  (
129  IOobject::groupName(alphaName_, U.group())
130  );
131 
132  return apply(AU, alpha, rho, U);
133  }
134 }
135 
136 
137 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual Type value(const scalar x) const =0
Return value as a function of (scalar) independent variable.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:482
static word groupName(Name name, const word &group)
List< label > labelList
A List of labels.
Definition: labelList.H:56
volScalarField scalarField(fieldObject, mesh)
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:76
static const Tensor I
Definition: Tensor.H:83
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:91
Field< tensor > tensorField
Specialisation of Field<T> for tensor.