A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise. More...
Classes | |
class | fvSolver |
Solver class returned by the solver function. More... | |
Public Member Functions | |
ClassName ("fvMatrix") | |
fvMatrix (const VolField< Type > &, const dimensionSet &) | |
Construct given a field to solve for. More... | |
fvMatrix (const fvMatrix< Type > &) | |
Copy constructor. More... | |
fvMatrix (const tmp< fvMatrix< Type >> &) | |
Copy constructor of tmp<fvMatrix<Type>> deleting argument. More... | |
fvMatrix (const VolField< Type > &, Istream &) | |
Construct from Istream given field to solve for. More... | |
tmp< fvMatrix< Type > > | clone () const |
Clone. More... | |
virtual | ~fvMatrix () |
Destructor. More... | |
VolField< Type > & | psi () |
const VolField< Type > & | psi () const |
const dimensionSet & | dimensions () const |
Field< Type > & | source () |
const Field< Type > & | source () const |
FieldField< Field, Type > & | internalCoeffs () |
fvBoundary scalar field containing pseudo-matrix coeffs More... | |
FieldField< Field, Type > & | boundaryCoeffs () |
fvBoundary scalar field containing pseudo-matrix coeffs More... | |
SurfaceField< Type > *& | faceFluxCorrectionPtr () |
Return pointer to face-flux non-orthogonal correction field. More... | |
template<template< class > class ListType> | |
void | setValues (const labelUList &cells, const ListType< Type > &values) |
Set solution in given cells to the specified values. More... | |
template<template< class > class ListType> | |
void | setValues (const labelUList &cells, const ListType< Type > &values, const scalarList &fractions, const bool hasDdt=true) |
Set solution in given cells to the specified values. More... | |
void | setReference (const label celli, const Type &value, const bool forceReference=false) |
Set reference level for solution. More... | |
void | setComponentReference (const label patchi, const label facei, const direction cmpt, const scalar value) |
Set reference level for a component of the solution. More... | |
scalar | relaxationFactor () const |
Return the relaxation factor for this equation. More... | |
void | relax (const scalar alpha) |
Relax matrix (for steady-state solution). More... | |
void | relax () |
Relax matrix (for steady-state solution). More... | |
void | boundaryManipulate (typename VolField< Type >::Boundary &values) |
Manipulate based on a boundary field. More... | |
autoPtr< fvSolver > | solver (const dictionary &) |
Construct and return the solver. More... | |
autoPtr< fvSolver > | solver () |
Construct and return the solver. More... | |
SolverPerformance< Type > | solve (const dictionary &) |
Solve segregated or coupled returning the solution statistics. More... | |
SolverPerformance< Type > | solveSegregated (const dictionary &) |
Solve segregated returning the solution statistics. More... | |
SolverPerformance< Type > | solveCoupled (const dictionary &) |
Solve coupled returning the solution statistics. More... | |
SolverPerformance< Type > | solve (const word &name) |
Solve segregated or coupled returning the solution statistics. More... | |
SolverPerformance< Type > | solve () |
Solve returning the solution statistics. More... | |
tmp< Field< Type > > | residual () const |
Return the matrix residual. More... | |
tmp< scalarField > | D () const |
Return the matrix scalar diagonal. More... | |
tmp< Field< Type > > | DD () const |
Return the matrix Type diagonal. More... | |
tmp< volScalarField > | A () const |
Return the central coefficient. More... | |
tmp< VolInternalField< Type > > | Su () const |
Return the explicit source. More... | |
tmp< volScalarField::Internal > | Sp () const |
Return the implicit source. More... | |
tmp< VolField< Type > > | H () const |
Return the H operation source. More... | |
tmp< volScalarField > | H1 () const |
Return H(1) More... | |
tmp< SurfaceField< Type > > | flux () const |
Return the face-flux field from the matrix. More... | |
void | operator= (const fvMatrix< Type > &) |
void | operator= (const tmp< fvMatrix< Type >> &) |
void | negate () |
void | operator+= (const fvMatrix< Type > &) |
void | operator+= (const tmp< fvMatrix< Type >> &) |
void | operator-= (const fvMatrix< Type > &) |
void | operator-= (const tmp< fvMatrix< Type >> &) |
void | operator+= (const DimensionedField< Type, volMesh > &) |
void | operator+= (const tmp< DimensionedField< Type, volMesh >> &) |
void | operator+= (const tmp< VolField< Type >> &) |
void | operator-= (const DimensionedField< Type, volMesh > &) |
void | operator-= (const tmp< DimensionedField< Type, volMesh >> &) |
void | operator-= (const tmp< VolField< Type >> &) |
void | operator+= (const dimensioned< Type > &) |
void | operator-= (const dimensioned< Type > &) |
void | operator+= (const zero &) |
void | operator-= (const zero &) |
void | operator*= (const volScalarField::Internal &) |
void | operator*= (const tmp< volScalarField::Internal > &) |
void | operator*= (const tmp< volScalarField > &) |
void | operator*= (const dimensioned< scalar > &) |
void | operator/= (const volScalarField::Internal &) |
void | operator/= (const tmp< volScalarField::Internal > &) |
void | operator/= (const tmp< volScalarField > &) |
void | operator/= (const dimensioned< scalar > &) |
void | setComponentReference (const label patchi, const label facei, const direction, const scalar value) |
Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolver > | solver (const dictionary &solverControls) |
Foam::solverPerformance | solveSegregated (const dictionary &solverControls) |
Foam::tmp< Foam::scalarField > | residual () const |
Foam::tmp< Foam::volScalarField > | H () const |
Foam::tmp< Foam::volScalarField > | H1 () const |
void | setComponentReference (const label patchi, const label facei, const direction, const scalar value) |
autoPtr< fvMatrix< scalar >::fvSolver > | solver (const dictionary &) |
solverPerformance | solveSegregated (const dictionary &) |
tmp< scalarField > | residual () const |
tmp< volScalarField > | H () const |
tmp< volScalarField > | H1 () const |
Public Member Functions inherited from refCount | |
int | count () const |
Return the current reference count. More... | |
bool | unique () const |
Return true if the reference count is zero. More... | |
void | operator++ () |
Increment the reference count. More... | |
void | operator++ (int) |
Increment the reference count. More... | |
void | operator-- () |
Decrement the reference count. More... | |
void | operator-- (int) |
Decrement the reference count. More... | |
void | operator= (const refCount &)=delete |
Disallow bitwise assignment. More... | |
Public Member Functions inherited from lduMatrix | |
ClassName ("lduMatrix") | |
lduMatrix (const lduMesh &) | |
Construct given an LDU addressed mesh. More... | |
lduMatrix (const lduMatrix &) | |
Copy constructor. More... | |
lduMatrix (lduMatrix &, bool reuse) | |
Copy constructor or reuse as specified. More... | |
lduMatrix (const lduMesh &, Istream &) | |
Construct given an LDU addressed mesh and an Istream. More... | |
~lduMatrix () | |
Destructor. More... | |
const lduMesh & | mesh () const |
Return the LDU mesh from which the addressing is obtained. More... | |
const lduAddressing & | lduAddr () const |
Return the LDU addressing. More... | |
const lduSchedule & | patchSchedule () const |
Return the patch evaluation schedule. More... | |
scalarField & | lower () |
scalarField & | diag () |
scalarField & | upper () |
scalarField & | lower (const label size) |
scalarField & | diag (const label nCoeffs) |
scalarField & | upper (const label nCoeffs) |
const scalarField & | lower () const |
const scalarField & | diag () const |
const scalarField & | upper () const |
bool | hasDiag () const |
bool | hasUpper () const |
bool | hasLower () const |
bool | diagonal () const |
bool | symmetric () const |
bool | asymmetric () const |
void | sumDiag () |
void | negSumDiag () |
void | sumMagOffDiag (scalarField &sumOff) const |
void | Amul (scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const |
Matrix multiplication with updated interfaces. More... | |
void | Tmul (scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const |
Matrix transpose multiplication with updated interfaces. More... | |
void | sumA (scalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const |
Sum the coefficients on each row of the matrix. More... | |
void | residual (scalarField &rA, const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const |
tmp< scalarField > | residual (const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const |
void | initMatrixInterfaces (const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const |
Initialise the update of interfaced interfaces. More... | |
void | updateMatrixInterfaces (const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const |
Update interfaced interfaces for matrix operations. More... | |
template<class Type > | |
tmp< Field< Type > > | H (const Field< Type > &) const |
template<class Type > | |
tmp< Field< Type > > | H (const tmp< Field< Type >> &) const |
tmp< scalarField > | H1 () const |
template<class Type > | |
tmp< Field< Type > > | faceH (const Field< Type > &) const |
template<class Type > | |
tmp< Field< Type > > | faceH (const tmp< Field< Type >> &) const |
InfoProxy< lduMatrix > | info () const |
Return info proxy. More... | |
void | operator= (const lduMatrix &) |
void | negate () |
void | operator+= (const lduMatrix &) |
void | operator-= (const lduMatrix &) |
void | operator*= (const scalarField &) |
void | operator*= (scalar) |
void | operator/= (const scalarField &) |
void | operator/= (scalar) |
template<class Type > | |
Foam::tmp< Foam::Field< Type > > | H (const Field< Type > &psi) const |
template<class Type > | |
Foam::tmp< Foam::Field< Type > > | H (const tmp< Field< Type >> &tpsi) const |
template<class Type > | |
Foam::tmp< Foam::Field< Type > > | faceH (const Field< Type > &psi) const |
template<class Type > | |
Foam::tmp< Foam::Field< Type > > | faceH (const tmp< Field< Type >> &tpsi) const |
Protected Member Functions | |
template<class Type2 > | |
void | addToInternalField (const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const |
Add patch contribution to internal field. More... | |
template<class Type2 > | |
void | addToInternalField (const labelUList &addr, const tmp< Field< Type2 >> &tpf, Field< Type2 > &intf) const |
template<class Type2 > | |
void | subtractFromInternalField (const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const |
Subtract patch contribution from internal field. More... | |
template<class Type2 > | |
void | subtractFromInternalField (const labelUList &addr, const tmp< Field< Type2 >> &tpf, Field< Type2 > &intf) const |
void | addBoundaryDiag (scalarField &diag, const direction cmpt) const |
void | addCmptAvBoundaryDiag (scalarField &diag) const |
void | addBoundarySource (Field< Type > &source, const bool couples=true) const |
void | setValue (const label celli, const Type &value) |
Set solution in the given cell to the specified value. More... | |
void | setValue (const label celli, const Type &value, const scalar fraction, const scalarField &ddtDiag) |
Set solution in the given cell to the specified value. More... | |
Protected Member Functions inherited from refCount | |
refCount () | |
Construct null initialising count to 0. More... | |
refCount (const refCount &)=delete | |
Disallow copy. More... | |
Friends | |
class | fvSolver |
Declare friendship with the fvSolver class. More... | |
tmp< VolField< Type > > | operator& (const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &) |
tmp< VolField< Type > > | operator& (const fvMatrix< Type > &, const tmp< VolField< Type >> &) |
tmp< VolField< Type > > | operator& (const tmp< fvMatrix< Type >> &, const DimensionedField< Type, volMesh > &) |
tmp< VolField< Type > > | operator& (const tmp< fvMatrix< Type >> &, const tmp< VolField< Type >> &) |
Ostream & | operator (Ostream &, const fvMatrix< Type > &) |
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition at line 114 of file fvMatrix.H.
fvMatrix | ( | const VolField< Type > & | psi, |
const dimensionSet & | ds | ||
) |
Construct given a field to solve for.
Definition at line 284 of file fvMatrix.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), Foam::endl(), regIOobject::eventNo(), forAll, InfoInFunction, patchi, fvMatrix< Type >::psi(), GeometricBoundaryField< Type, PatchField, GeoMesh >::updateCoeffs(), and Foam::Zero.
Copy constructor.
Definition at line 339 of file fvMatrix.C.
References Foam::endl(), and InfoInFunction.
Copy constructor of tmp<fvMatrix<Type>> deleting argument.
Definition at line 368 of file fvMatrix.C.
References Foam::endl(), and InfoInFunction.
Construct from Istream given field to solve for.
Definition at line 422 of file fvMatrix.C.
References Foam::endl(), forAll, InfoInFunction, patchi, fvMatrix< Type >::psi(), and Foam::Zero.
|
virtual |
Destructor.
Definition at line 482 of file fvMatrix.C.
References Foam::endl(), and InfoInFunction.
|
protected |
Add patch contribution to internal field.
Definition at line 39 of file fvMatrix.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, UList< T >::size(), and List< T >::size().
|
protected |
Definition at line 62 of file fvMatrix.C.
|
protected |
Subtract patch contribution from internal field.
Definition at line 76 of file fvMatrix.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, UList< T >::size(), and List< T >::size().
|
protected |
Definition at line 99 of file fvMatrix.C.
|
protected |
Definition at line 112 of file fvMatrix.C.
References Foam::component(), Foam::diag(), forAll, and patchi.
Referenced by fvMatrix< Type >::residual(), and fvMatrix< Type >::solveSegregated().
|
protected |
Definition at line 131 of file fvMatrix.C.
References Foam::cmptAv(), Foam::diag(), forAll, and patchi.
|
protected |
Definition at line 146 of file fvMatrix.C.
References Foam::cmptMultiply(), fvPatchField< Type >::coupled(), forAll, patchi, and fvPatchField< Type >::patchNeighbourField().
Referenced by fvMatrix< Type >::H(), fvMatrix< Type >::residual(), and fvMatrix< Type >::solveSegregated().
|
protected |
Set solution in the given cell to the specified value.
Definition at line 178 of file fvMatrix.C.
References Foam::constant::universal::c, cells, primitiveMesh::cells(), Foam::diag(), forAll, primitiveMesh::isInternalFace(), fvMesh::mesh(), fvMesh::neighbour(), primitiveMesh::nInternalFaces(), fvMesh::owner(), patches, fvMesh::polyBFacePatches(), fvMesh::polyBFacePatchFaces(), primitiveFieldRef(), psi, and Foam::Zero.
|
protected |
Set solution in the given cell to the specified value.
Definition at line 258 of file fvMatrix.C.
References Foam::diag(), primitiveFieldRef(), and psi.
ClassName | ( | "fvMatrix< Type >" | ) |
Foam::tmp< Foam::fvMatrix< Type > > clone |
Clone.
Definition at line 470 of file fvMatrix.C.
|
inline |
Definition at line 289 of file fvMatrix.H.
Referenced by fvTotalSource::addSource(), interfaceTurbulenceDamping::addSup(), VoFSolidificationMelting::addSup(), VoFCavitation::addSup(), VoFClouds::addSup(), VoFFilmTransfer::addSup(), homogeneousCondensation::addSup(), homogeneousLiquidPhaseSeparation::addSup(), coefficientPhaseChange::addSup(), coefficientMassTransfer::addSup(), waveForcing::addSup(), filmCloudTransfer::addSup(), filmVoFTransfer::addSup(), VoFTurbulenceDamping::addSup(), solidThermalEquilibrium::addSup(), radiation::addSup(), clouds::addSup(), porosityForce::addSup(), solidificationMelting::addSup(), interRegionPorosityForce::addSup(), massSourceBase::addSup(), massTransfer::addSup(), volumeSource::addSup(), fvSpecificSource::addSupType(), fvTotalSource::addSupType(), massTransfer::addSupType(), Foam::checkMethod(), fvConstraints::constrain(), fvMatrix< Type >::fvMatrix(), zeroDimensionalFixedPressureConstraint::pEqnSource(), fvMatrix< Type >::fvSolver::solve(), and fvMatrix< Type >::solveSegregated().
|
inline |
Definition at line 297 of file fvMatrix.H.
|
inline |
Definition at line 302 of file fvMatrix.H.
Referenced by fvTotalSource::addSource(), porosityForce::addSup(), interRegionPorosityForce::addSup(), rotorDisk::addSup(), fvTotalSource::addSupType(), massTransfer::addSupType(), Foam::checkMethod(), zeroDimensionalFixedPressureConstraint::constrain(), meanVelocityForce::constrain(), fvMatrix< Type >::operator*=(), fvMatrix< Type >::operator/=(), and zeroDimensionalFixedPressureConstraint::pEqnSource().
|
inline |
Definition at line 307 of file fvMatrix.H.
Referenced by fvTotalSource::addSource(), propellerDisk::addSup(), heatSource::addSup(), effectivenessHeatExchanger::addSup(), actuationDisk::addSup(), radialActuationDisk::addSup(), solidificationMelting::addSup(), fvTotalSource::addSupType(), EulerD2dt2Scheme< Type >::fvmD2dt2(), backwardDdtScheme< Type >::fvmDdt(), CoEulerDdtScheme< Type >::fvmDdt(), CrankNicolsonDdtScheme< Type >::fvmDdt(), EulerDdtScheme< Type >::fvmDdt(), localEulerDdtScheme< Type >::fvmDdt(), SLTSDdtScheme< Type >::fvmDdt(), gaussLaplacianScheme< Type, GType >::fvmLaplacian(), isothermalFilm::momentumPredictor(), incompressibleDenseParticleFluid::prePredictor(), ThermoCloud< CloudType >::Sh(), ReactingCloud< CloudType >::Srho(), MomentumCloud< CloudType >::SU(), and ReactingCloud< CloudType >::SYi().
|
inline |
Definition at line 312 of file fvMatrix.H.
|
inline |
fvBoundary scalar field containing pseudo-matrix coeffs
for internal cells
Definition at line 319 of file fvMatrix.H.
Referenced by gaussConvectionScheme< Type >::fvmDiv(), and gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected().
|
inline |
fvBoundary scalar field containing pseudo-matrix coeffs
for boundary cells
Definition at line 326 of file fvMatrix.H.
Referenced by gaussConvectionScheme< Type >::fvmDiv(), and gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected().
|
inline |
Return pointer to face-flux non-orthogonal correction field.
Definition at line 332 of file fvMatrix.H.
Referenced by gaussConvectionScheme< Type >::fvmDiv(), and gaussLaplacianScheme< Type, GType >::fvmLaplacian().
void setValues | ( | const labelUList & | cells, |
const ListType< Type > & | values | ||
) |
Set solution in given cells to the specified values.
and eliminate the corresponding equations from the matrix.
Definition at line 501 of file fvMatrix.C.
References forAll.
Referenced by fixedTemperature::constrain(), and fixedInternalValueFvPatchField< Type >::manipulateMatrix().
void setValues | ( | const labelUList & | cells, |
const ListType< Type > & | values, | ||
const scalarList & | fractions, | ||
const bool | hasDdt = true |
||
) |
Set solution in given cells to the specified values.
and eliminate the corresponding equations from the matrix.
Definition at line 517 of file fvMatrix.C.
References alpha(), Foam::fvm::ddt(), lduMatrix::diag(), Foam::diag(), forAll, fractions(), and lduMatrix::hasDiag().
void setReference | ( | const label | celli, |
const Type & | value, | ||
const bool | forceReference = false |
||
) |
Set reference level for solution.
Definition at line 562 of file fvMatrix.C.
References Foam::diag().
Referenced by Foam::fv::CorrectPhi(), incompressibleDenseParticleFluid::correctPressure(), incompressibleFluid::correctPressure(), twoPhaseSolver::incompressiblePressureCorrector(), main(), and incompressibleMultiphaseVoF::pressureCorrector().
void setComponentReference | ( | const label | patchi, |
const label | facei, | ||
const direction | cmpt, | ||
const scalar | value | ||
) |
Set reference level for a component of the solution.
on a given patch face
Definition at line 33 of file fvMatrixSolve.C.
References Foam::diag(), and patchi.
Foam::scalar relaxationFactor |
Return the relaxation factor for this equation.
Definition at line 578 of file fvMatrix.C.
References solutionControl::finalIteration().
void relax | ( | const scalar | alpha | ) |
Relax matrix (for steady-state solution).
alpha = 1 : diagonally equal alpha < 1 : diagonally dominant alpha = 0 : do nothing Note: Requires positive diagonal.
Definition at line 603 of file fvMatrix.C.
References alpha(), Foam::cmptMag(), Foam::cmptMax(), Foam::cmptMin(), Foam::component(), fvPatchField< Type >::coupled(), D, Foam::diag(), Foam::endl(), forAll, InfoInFunction, Foam::mag(), Foam::max(), UPstream::msgType(), Foam::nl, patchi, Foam::reduce(), Foam::returnReduce(), Foam::fvm::S(), and List< T >::size().
Referenced by XiFluid::bSolve(), kineticTheoryModel::correct(), IATE::correct(), fractal::correct(), Maxwell< BasicMomentumTransportModel >::correct(), dynamicLagrangian< BasicMomentumTransportModel >::correct(), radiativeIntensityRay::correct(), advectionDiffusion::correct(), age::execute(), phaseScalarTransport::execute(), scalarTransport::execute(), main(), incompressibleDenseParticleFluid::momentumPredictor(), populationBalanceModel::solve(), compressibleMultiphaseVoF::thermophysicalPredictor(), compressibleVoF::thermophysicalPredictor(), film::thermophysicalPredictor(), multicomponentFluid::thermophysicalPredictor(), and solid::thermophysicalPredictor().
void relax |
Relax matrix (for steady-state solution).
alpha is read from controlDict
Definition at line 754 of file fvMatrix.C.
References relax().
void boundaryManipulate | ( | typename VolField< Type >::Boundary & | values | ) |
Manipulate based on a boundary field.
Definition at line 761 of file fvMatrix.C.
autoPtr<fvSolver> solver | ( | const dictionary & | ) |
Construct and return the solver.
Use the given solver controls
Foam::autoPtr< typename Foam::fvMatrix< Type >::fvSolver > solver |
Construct and return the solver.
Solver controls read from fvSolution
Definition at line 281 of file fvMatrixSolve.C.
Foam::SolverPerformance< Type > solve | ( | const dictionary & | solverControls | ) |
Solve segregated or coupled returning the solution statistics.
Use the given solver controls
Definition at line 57 of file fvMatrixSolve.C.
References Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, Foam::Info, dictionary::lookupOrDefault(), dictionary::readIfPresent(), and Foam::type().
Referenced by XiFluid::bSolve(), Implicit< CloudType >::cacheFields(), kineticTheoryModel::correct(), IATE::correct(), fractal::correct(), Maxwell< BasicMomentumTransportModel >::correct(), dynamicLagrangian< BasicMomentumTransportModel >::correct(), advectionDiffusion::correct(), Foam::fv::CorrectPhi(), incompressibleDenseParticleFluid::correctPressure(), incompressibleFluid::correctPressure(), XiFluid::EauSolve(), age::execute(), phaseScalarTransport::execute(), scalarTransport::execute(), XiFluid::ftSolve(), Foam::hydrostaticInitialisation(), twoPhaseSolver::incompressiblePressureCorrector(), main(), incompressibleDriftFlux::prePredictor(), compressibleMultiphaseVoF::pressureCorrector(), compressibleVoF::pressureCorrector(), incompressibleMultiphaseVoF::pressureCorrector(), solidDisplacement::pressureCorrector(), populationBalanceModel::solve(), phaseSystem::solve(), compressibleMultiphaseVoF::thermophysicalPredictor(), compressibleVoF::thermophysicalPredictor(), film::thermophysicalPredictor(), multicomponentFluid::thermophysicalPredictor(), and solid::thermophysicalPredictor().
Foam::SolverPerformance< Type > solveSegregated | ( | const dictionary & | solverControls | ) |
Solve segregated returning the solution statistics.
Use the given solver controls
Definition at line 104 of file fvMatrixSolve.C.
References Field< Type >::component(), Foam::diag(), Foam::endl(), Foam::Info, Foam::compressible::New(), SolverPerformance< Type >::print(), psi, GeometricField< Type, PatchField, GeoMesh >::replace(), and SolverPerformance< Type >::solverName().
Foam::SolverPerformance< Type > solveCoupled | ( | const dictionary & | solverControls | ) |
Solve coupled returning the solution statistics.
Use the given solver controls
Definition at line 220 of file fvMatrixSolve.C.
References GeometricField< Type, PatchField, GeoMesh >::component(), Foam::diag(), Foam::endl(), Foam::Info, SolverPerformance< Type >::print(), and psi.
Foam::SolverPerformance< Type > solve | ( | const word & | name | ) |
Solve segregated or coupled returning the solution statistics.
Solver controls read from fvSolution
Definition at line 315 of file fvMatrixSolve.C.
References solutionControl::finalIteration(), Foam::name(), and fvMatrix< Type >::solve().
Foam::SolverPerformance< Type > solve |
Solve returning the solution statistics.
Solver controls read from fvSolution
Definition at line 331 of file fvMatrixSolve.C.
References solutionControl::finalIteration().
Referenced by fvMatrix< Type >::solve().
Foam::tmp< Foam::Field< Type > > residual |
Return the matrix residual.
Definition at line 348 of file fvMatrixSolve.C.
References fvMatrix< Type >::addBoundaryDiag(), fvMatrix< Type >::addBoundarySource(), Field< Type >::component(), tmp< T >::ref(), Field< Type >::replace(), and lduMatrix::residual().
Return the matrix scalar diagonal.
Definition at line 775 of file fvMatrix.C.
References Foam::diag(), and tmp< T >::ref().
Foam::tmp< Foam::Field< Type > > DD |
Return the matrix Type diagonal.
Definition at line 784 of file fvMatrix.C.
References fvPatchField< Type >::coupled(), Foam::diag(), forAll, patchi, tmp< T >::ref(), and List< T >::size().
Return the central coefficient.
Definition at line 808 of file fvMatrix.C.
References D, Foam::dimVolume, GeometricField< Type, PatchField, GeoMesh >::New(), and tmp< T >::ref().
Referenced by meanVelocityForce::constrain(), incompressibleDenseParticleFluid::correctPressure(), MomentumTransferPhaseSystem< BasePhaseSystem >::invADVfs(), MomentumTransferPhaseSystem< BasePhaseSystem >::invADVs(), main(), and compressibleVoF::thermophysicalPredictor().
Foam::tmp< Foam::VolInternalField< Type > > Su |
Return the explicit source.
Definition at line 829 of file fvMatrix.C.
References Foam::dimVolume.
Referenced by compressibleVoF::alphaSuSp(), incompressibleDriftFlux::alphaSuSp(), and incompressibleVoF::alphaSuSp().
Return the implicit source.
Definition at line 847 of file fvMatrix.C.
References Foam::diag(), Foam::dimVolume, and DimensionedField< Type, GeoMesh >::New().
Referenced by compressibleVoF::alphaSuSp(), incompressibleDriftFlux::alphaSuSp(), and incompressibleVoF::alphaSuSp().
Foam::tmp< Foam::VolField< Type > > H |
Return the H operation source.
Definition at line 868 of file fvMatrix.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::dimVolume, lduMatrix::H(), Field< Type >::negate(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), tmp< T >::ref(), GeometricField< Type, PatchField, GeoMesh >::replace(), and Field< Type >::replace().
Referenced by incompressibleDenseParticleFluid::correctPressure(), MomentumTransferPhaseSystem< BasePhaseSystem >::invADVfs(), MomentumTransferPhaseSystem< BasePhaseSystem >::invADVs(), and main().
Return H(1)
Definition at line 923 of file fvMatrix.C.
References Foam::component(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), fvPatchField< Type >::coupled(), Foam::dimVolume, forAll, lduMatrix::H1(), GeometricField< Type, PatchField, GeoMesh >::New(), patchi, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), tmp< T >::ref(), and List< T >::size().
Foam::tmp< Foam::SurfaceField< Type > > flux |
Return the face-flux field from the matrix.
Definition at line 963 of file fvMatrix.C.
References Foam::abort(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), Foam::cmptMultiply(), Foam::dimensions(), lduMatrix::faceH(), Foam::FatalError, FatalErrorInFunction, forAll, patchi, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), tmp< T >::ref(), and Field< Type >::replace().
Referenced by Implicit< CloudType >::cacheFields(), Foam::fv::CorrectPhi(), incompressibleDenseParticleFluid::correctPressure(), incompressibleFluid::correctPressure(), twoPhaseSolver::incompressiblePressureCorrector(), main(), incompressibleDriftFlux::prePredictor(), compressibleMultiphaseVoF::pressureCorrector(), compressibleVoF::pressureCorrector(), incompressibleMultiphaseVoF::pressureCorrector(), solidDisplacement::pressureCorrector(), and phaseSystem::solve().
Definition at line 1044 of file fvMatrix.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and lduMatrix::operator=().
Definition at line 1080 of file fvMatrix.C.
void negate |
Definition at line 1088 of file fvMatrix.C.
References lduMatrix::negate().
Definition at line 1103 of file fvMatrix.C.
References Foam::checkMethod(), and lduMatrix::operator+=().
Definition at line 1129 of file fvMatrix.C.
Definition at line 1137 of file fvMatrix.C.
References Foam::checkMethod(), and lduMatrix::operator-=().
Definition at line 1161 of file fvMatrix.C.
void operator+= | ( | const DimensionedField< Type, volMesh > & | su | ) |
Definition at line 1169 of file fvMatrix.C.
References Foam::checkMethod().
void operator+= | ( | const tmp< DimensionedField< Type, volMesh >> & | tsu | ) |
Definition at line 1180 of file fvMatrix.C.
Definition at line 1191 of file fvMatrix.C.
void operator-= | ( | const DimensionedField< Type, volMesh > & | su | ) |
Definition at line 1202 of file fvMatrix.C.
References Foam::checkMethod().
void operator-= | ( | const tmp< DimensionedField< Type, volMesh >> & | tsu | ) |
Definition at line 1213 of file fvMatrix.C.
Definition at line 1224 of file fvMatrix.C.
void operator+= | ( | const dimensioned< Type > & | su | ) |
Definition at line 1235 of file fvMatrix.C.
References DimensionedField< Type, GeoMesh >::mesh(), and psi.
void operator-= | ( | const dimensioned< Type > & | su | ) |
Definition at line 1245 of file fvMatrix.C.
References DimensionedField< Type, GeoMesh >::mesh(), and psi.
Definition at line 1255 of file fvMatrix.C.
Definition at line 1263 of file fvMatrix.C.
void operator*= | ( | const volScalarField::Internal & | dsf | ) |
Definition at line 1271 of file fvMatrix.C.
References Foam::abort(), fvMatrix< Type >::dimensions(), Foam::FatalError, FatalErrorInFunction, forAll, lduMatrix::operator*=(), and patchi.
void operator*= | ( | const tmp< volScalarField::Internal > & | tdsf | ) |
Definition at line 1304 of file fvMatrix.C.
void operator*= | ( | const tmp< volScalarField > & | tvsf | ) |
Definition at line 1315 of file fvMatrix.C.
void operator*= | ( | const dimensioned< scalar > & | ds | ) |
Definition at line 1326 of file fvMatrix.C.
References fvMatrix< Type >::dimensions(), and lduMatrix::operator*=().
void operator/= | ( | const volScalarField::Internal & | dsf | ) |
Definition at line 1345 of file fvMatrix.C.
References Foam::abort(), fvMatrix< Type >::dimensions(), Foam::FatalError, FatalErrorInFunction, forAll, lduMatrix::operator/=(), and patchi.
void operator/= | ( | const tmp< volScalarField::Internal > & | tdsf | ) |
Definition at line 1378 of file fvMatrix.C.
void operator/= | ( | const tmp< volScalarField > & | tvsf | ) |
Definition at line 1389 of file fvMatrix.C.
void operator/= | ( | const dimensioned< scalar > & | ds | ) |
Definition at line 1400 of file fvMatrix.C.
References fvMatrix< Type >::dimensions(), and lduMatrix::operator/=().
void setComponentReference | ( | const label | patchi, |
const label | facei, | ||
const | direction, | ||
const scalar | value | ||
) |
Definition at line 33 of file fvScalarMatrix.C.
References Foam::diag(), and patchi.
Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolver > solver | ( | const dictionary & | solverControls | ) |
Definition at line 58 of file fvScalarMatrix.C.
References Foam::diag(), Foam::endl(), Foam::Info, and Foam::compressible::New().
Foam::solverPerformance solveSegregated | ( | const dictionary & | solverControls | ) |
Definition at line 138 of file fvScalarMatrix.C.
References fvMatrix< Type >::addBoundaryDiag(), fvMatrix< Type >::addBoundarySource(), Residuals< Type >::append(), lduMatrix::diag(), Foam::endl(), Foam::Info, lduMatrix::mesh(), lduMatrix::solver::New(), SolverPerformance< Type >::print(), and fvMatrix< Type >::psi().
Foam::tmp< Foam::scalarField > residual | ( | ) | const |
Definition at line 188 of file fvScalarMatrix.C.
References fvMatrix< Type >::addBoundaryDiag(), fvMatrix< Type >::addBoundarySource(), tmp< T >::ref(), and lduMatrix::residual().
Foam::tmp< Foam::volScalarField > H | ( | ) | const |
Definition at line 212 of file fvScalarMatrix.C.
References fvMatrix< Type >::addBoundarySource(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimVolume, lduMatrix::H(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and tmp< T >::ref().
Foam::tmp< Foam::volScalarField > H1 | ( | ) | const |
Definition at line 237 of file fvScalarMatrix.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimVolume, lduMatrix::H1(), GeometricField< Type, PatchField, GeoMesh >::New(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and tmp< T >::ref().
void setComponentReference | ( | const label | patchi, |
const label | facei, | ||
const | direction, | ||
const scalar | value | ||
) |
autoPtr< fvMatrix< scalar >::fvSolver > solver | ( | const dictionary & | ) |
solverPerformance solveSegregated | ( | const dictionary & | ) |
tmp< scalarField > residual | ( | ) | const |
tmp< volScalarField > H | ( | ) | const |
tmp< volScalarField > H1 | ( | ) | const |
|
friend |
Declare friendship with the fvSolver class.
Definition at line 147 of file fvMatrix.H.
|
friend |
|
friend |
|
friend |