solidificationMeltingSourceTemplates.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) 2014-2021 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 "fvMatrices.H"
27 #include "fvcDdt.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class RhoFieldType>
32 void Foam::fv::solidificationMeltingSource::apply
33 (
34  const RhoFieldType& rho,
35  fvMatrix<scalar>& eqn
36 ) const
37 {
38  if (debug)
39  {
40  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
41  }
42 
43  const volScalarField Cp(this->Cp());
44 
45  update(Cp);
46 
48 
49  // Contributions added to rhs of solver equation
50  if (eqn.psi().dimensions() == dimTemperature)
51  {
52  eqn -= L/Cp*(fvc::ddt(rho, alpha1_));
53  }
54  else
55  {
56  eqn -= L*(fvc::ddt(rho, alpha1_));
57  }
58 }
59 
60 
61 // ************************************************************************* //
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the first temporal derivative.
const dimensionSet dimEnergy
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
const dimensionSet dimTemperature