faMatrixSolve.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-2017 Wikki Ltd
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Description
28  Finite-Area matrix basic solvers.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "PrecisionAdaptor.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const label patchi,
40  const label facei,
41  const direction cmpt,
42  const scalar value
43 )
44 {
45  internalCoeffs_[patchi][facei].component(cmpt) +=
46  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
47 
48  boundaryCoeffs_[patchi][facei].component(cmpt) +=
49  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]*value;
50 }
51 
52 
53 template<class Type>
55 (
56  const dictionary& solverControls
57 )
58 {
60  << "solving faMatrix<Type>"
61  << endl;
62 
63  SolverPerformance<Type> solverPerfVec
64  (
65  "faMatrix<Type>::solve",
66  psi_.name()
67  );
68 
69  scalarField saveDiag(diag());
70 
71  Field<Type> source(source_);
72  addBoundarySource(source);
73 
74  // Note: make a copy of interfaces: no longer a reference
75  lduInterfaceFieldPtrsList interfaces =
76  psi_.boundaryField().scalarInterfaces();
77 
78  // Cast into a non-const to solve
81 
82  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
83  {
84  // copy field and source
85 
86  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
87  addBoundaryDiag(diag(), cmpt);
88 
89  scalarField sourceCmpt(source.component(cmpt));
90 
91  FieldField<Field, scalar> bouCoeffsCmpt
92  (
93  boundaryCoeffs_.component(cmpt)
94  );
95 
96  FieldField<Field, scalar> intCoeffsCmpt
97  (
98  internalCoeffs_.component(cmpt)
99  );
100 
101  // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
102  // bouCoeffsCmpt for the explicit part of the coupled boundary
103  // conditions
104  {
105  PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
106  ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
107 
108  initMatrixInterfaces
109  (
110  true,
111  bouCoeffsCmpt,
112  interfaces,
113  psiCmpt_ss(),
114  sourceCmpt_ss.ref(),
115  cmpt
116  );
117 
118  updateMatrixInterfaces
119  (
120  true,
121  bouCoeffsCmpt,
122  interfaces,
123  psiCmpt_ss(),
124  sourceCmpt_ss.ref(),
125  cmpt
126  );
127  }
128 
129  solverPerformance solverPerf;
130 
131  // Solver call
132  solverPerf = lduMatrix::solver::New
133  (
134  psi_.name() + pTraits<Type>::componentNames[cmpt],
135  *this,
136  bouCoeffsCmpt,
137  intCoeffsCmpt,
138  interfaces,
139  solverControls
140  )->solve(psiCmpt, sourceCmpt, cmpt);
141 
143  {
144  solverPerf.print(Info);
145  }
146 
147  solverPerfVec.replace(cmpt, solverPerf);
148  solverPerfVec.solverName() = solverPerf.solverName();
149 
150  psi.primitiveFieldRef().replace(cmpt, psiCmpt);
151  diag() = saveDiag;
152  }
153 
154  psi.correctBoundaryConditions();
155 
156  psi.mesh().setSolverPerformance(psi.name(), solverPerfVec);
157 
158  return solverPerfVec;
159 }
160 
161 
162 template<class Type>
164 {
165  return solve(faMat_.psi().mesh().solverDict(faMat_.psi().name()));
166 }
167 
168 
169 template<class Type>
171 {
172  return solve(this->psi().mesh().solverDict(this->psi().name()));
173 }
174 
175 
176 template<class Type>
178 {
179  tmp<Field<Type>> tres(source_);
180  Field<Type>& res = tres().ref();
181 
182  addBoundarySource(res);
183 
184  lduInterfaceFieldPtrsList interfaces =
185  psi_.boundaryField().scalarInterfaces();
186 
187  // Loop over field components
188  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
189  {
190  scalarField psiCmpt(psi_.internalField().component(cmpt));
191 
192  scalarField boundaryDiagCmpt(psi_.size(), Zero);
193  addBoundaryDiag(boundaryDiagCmpt, cmpt);
194 
195  FieldField<Field, scalar> bouCoeffsCmpt
196  (
197  boundaryCoeffs_.component(cmpt)
198  );
199 
200  res.replace
201  (
202  cmpt,
204  (
205  psiCmpt,
206  res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
207  bouCoeffsCmpt,
208  interfaces,
209  cmpt
210  )
211  );
212  }
213 
214  return tres;
215 }
216 
217 
218 // ************************************************************************* //
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
PrecisionAdaptor.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::faMatrix::setComponentReference
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution.
Definition: faMatrixSolve.C:38
Foam::SolverPerformance::solverName
const word & solverName() const
Return solver name.
Definition: SolverPerformance.H:149
solve
CEqn solve()
Foam::SolverPerformance::print
void print(Ostream &os) const
Print summary of solver performance to the given stream.
Definition: SolverPerformance.C:99
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:356
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::UPtrList< const lduInterfaceField >
Foam::faMatrix::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: faMatrixSolve.C:170
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:574
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:562
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::PrecisionAdaptor
A Field wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:158
Foam::PrecisionAdaptor::ref
FieldType & ref()
Allow modification without const-ref check.
Definition: PrecisionAdaptor.H:219
Foam::faMatrix::psi
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:221
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::ConstPrecisionAdaptor
A const Field wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:50
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::faMatrix::residual
tmp< Field< Type > > residual() const
Return the matrix residual.
Definition: faMatrixSolve.C:177
Foam::faMatrix::faSolver::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: faMatrixSolve.C:163
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::lduMatrix::residual
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
Definition: lduMatrixATmul.C:211
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51