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  auto& psi =
65 
66  SolverPerformance<Type> solverPerfVec
67  (
68  "faMatrix<Type>::solve",
69  psi.name()
70  );
71 
72  scalarField saveDiag(diag());
73 
74  Field<Type> source(source_);
75  addBoundarySource(source);
76 
77  // Note: make a copy of interfaces: no longer a reference
78  lduInterfaceFieldPtrsList interfaces =
79  psi_.boundaryField().scalarInterfaces();
80 
81  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
82  {
83  // copy field and source
84 
85  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
86  addBoundaryDiag(diag(), cmpt);
87 
88  scalarField sourceCmpt(source.component(cmpt));
89 
90  FieldField<Field, scalar> bouCoeffsCmpt
91  (
92  boundaryCoeffs_.component(cmpt)
93  );
94 
95  FieldField<Field, scalar> intCoeffsCmpt
96  (
97  internalCoeffs_.component(cmpt)
98  );
99 
100  // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
101  // bouCoeffsCmpt for the explicit part of the coupled boundary
102  // conditions
103  {
104  PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
105  ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
106 
107  const label startRequest = Pstream::nRequests();
108 
109  initMatrixInterfaces
110  (
111  true,
112  bouCoeffsCmpt,
113  interfaces,
114  psiCmpt_ss(),
115  sourceCmpt_ss.ref(),
116  cmpt
117  );
118 
119  updateMatrixInterfaces
120  (
121  true,
122  bouCoeffsCmpt,
123  interfaces,
124  psiCmpt_ss(),
125  sourceCmpt_ss.ref(),
126  cmpt,
127  startRequest
128  );
129  }
130 
131  solverPerformance solverPerf;
132 
133  // Solver call
134  solverPerf = lduMatrix::solver::New
135  (
136  psi_.name() + pTraits<Type>::componentNames[cmpt],
137  *this,
138  bouCoeffsCmpt,
139  intCoeffsCmpt,
140  interfaces,
141  solverControls
142  )->solve(psiCmpt, sourceCmpt, cmpt);
143 
145  {
146  solverPerf.print(Info);
147  }
148 
149  solverPerfVec.replace(cmpt, solverPerf);
150  solverPerfVec.solverName() = solverPerf.solverName();
151 
152  psi.primitiveFieldRef().replace(cmpt, psiCmpt);
153  diag() = saveDiag;
154  }
155 
156  psi.correctBoundaryConditions();
157 
158  psi.mesh().setSolverPerformance(psi.name(), solverPerfVec);
159 
160  return solverPerfVec;
161 }
162 
163 
164 template<class Type>
166 {
167  return solve(faMat_.psi().mesh().solverDict(faMat_.psi().name()));
168 }
169 
170 
171 template<class Type>
173 {
174  return solve(this->psi().mesh().solverDict(this->psi().name()));
175 }
176 
177 
178 template<class Type>
180 {
181  tmp<Field<Type>> tres(new Field<Type>(source_));
182  Field<Type>& res = tres().ref();
183 
184  addBoundarySource(res);
185 
186  // Loop over field components
187  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
188  {
189  scalarField psiCmpt(psi_.internalField().component(cmpt));
190 
191  scalarField boundaryDiagCmpt(psi_.size(), Zero);
192  addBoundaryDiag(boundaryDiagCmpt, cmpt);
193 
194  FieldField<Field, scalar> bouCoeffsCmpt
195  (
196  boundaryCoeffs_.component(cmpt)
197  );
198 
199  res.replace
200  (
201  cmpt,
203  (
204  psiCmpt,
205  res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
206  bouCoeffsCmpt,
207  psi_.boundaryField().scalarInterfaces(),
208  cmpt
209  )
210  );
211  }
212 
213  return tres;
214 }
215 
216 
217 // ************************************************************************* //
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:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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:350
Foam::faMatrix::setComponentReference
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
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::faMatrix::addBoundaryDiag
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: faMatrix.C:117
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
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:365
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:172
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:563
Foam::faMatrix::addBoundarySource
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: faMatrix.C:151
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:551
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 non-const Field/List wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:160
Foam::faMatrix::psi
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:239
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
Foam::ConstPrecisionAdaptor
A const Field/List wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:51
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::faMatrix::residual
tmp< Field< Type > > residual() const
Return the matrix residual.
Definition: faMatrixSolve.C:179
Foam::refPtr< Container< Type > >::ref
Container< Type > & ref() const
Definition: refPtrI.H:204
Foam::faMatrix::faSolver::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: faMatrixSolve.C:165
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:217
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51