Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
Backofen, Rainer
amdis
Commits
48245c1c
Commit
48245c1c
authored
Nov 10, 2012
by
Thomas Witkowski
Browse files
Fixed memory issues in navier stokes solver.
parent
8f6cf14b
Changes
6
Hide whitespace changes
Inline
Side-by-side
AMDiS/src/RefinementManager3d.cc
View file @
48245c1c
...
...
@@ -657,12 +657,12 @@ namespace AMDiS {
#ifdef HAVE_PARALLEL_DOMAIN_AMDIS
vector
<
FixRefinementPatch
::
EdgeInEl
>
refineEdges
;
FixRefinementPatch
::
getOtherEl
(
stack
,
refineEdges
);
if
(
refineEdges
.
size
())
{
MSG
(
"FIX REFINEMENT PATH ON ELEMENT %d %d: %d additional edges
\n
"
,
elInfo
->
getElement
()
->
getIndex
(),
elInfo
->
getMacroElement
()
->
getIndex
(),
refineEdges
.
size
());
}
//
if (refineEdges.size()) {
//
MSG("FIX REFINEMENT PATH ON ELEMENT %d %d: %d additional edges\n",
//
elInfo->getElement()->getIndex(),
//
elInfo->getMacroElement()->getIndex(),
//
refineEdges.size());
//
}
#endif
// === Traverse and refine the refinement patch. ====
...
...
@@ -677,15 +677,8 @@ namespace AMDiS {
// === If the refinement edge must be fixed, add also the other part ===
// === of this edge to the refinement patch. ===
for
(
int
edgeIndex
=
0
;
edgeIndex
<
static_cast
<
unsigned
int
>
(
refineEdges
.
size
());
edgeIndex
++
)
{
MSG
(
" IN REF FIX: %d %d
\n
"
,
refineEdges
[
edgeIndex
].
first
->
getIndex
(),
refineEdges
[
edgeIndex
].
second
);
}
for
(
int
edgeIndex
=
0
;
edgeIndex
<
static_cast
<
unsigned
int
>
(
refineEdges
.
size
());
edgeIndex
++
)
{
MSG
(
" MAKE -> %d
\n
"
,
edgeIndex
);
Element
*
otherEl
=
refineEdges
[
edgeIndex
].
first
;
TraverseStack
stack2
;
ElInfo
*
elInfo2
=
...
...
AMDiS/src/parallel/MeshDistributor.cc
View file @
48245c1c
...
...
@@ -1049,7 +1049,6 @@ namespace AMDiS {
// === Check the boundaries and adapt mesh if necessary. ===
MSG_DBG
(
"Run checkAndAdaptBoundary ...
\n
"
);
double
second
=
MPI
::
Wtime
();
// Check for periodic boundaries within rank's subdomain.
for
(
InteriorBoundary
::
iterator
it
(
intBoundary
.
getPeriodic
());
!
it
.
end
();
++
it
)
{
...
...
@@ -1064,12 +1063,9 @@ namespace AMDiS {
}
}
}
MSG
(
" -> create mesh structure codes needed %.5f seconds
\n
"
,
MPI
::
Wtime
()
-
second
);
second
=
MPI
::
Wtime
();
meshChanged
|=
checkAndAdaptBoundary
(
allBound
);
MSG
(
" -> checkAndAdaptBoundary needed %.5f seconds
\n
"
,
MPI
::
Wtime
()
-
second
);
// === Check on all ranks if at least one rank's mesh has changed. ===
...
...
AMDiS/src/parallel/PetscHelper.cc
View file @
48245c1c
...
...
@@ -249,13 +249,14 @@ namespace AMDiS {
}
void
setSolver
(
KSP
ksp
,
const
char
*
kspPrefix
,
KSPType
kspType
,
PCType
pcType
,
PetscReal
rtol
,
PetscReal
atol
,
PetscInt
maxIt
)
void
setSolverWithLu
(
KSP
ksp
,
const
char
*
kspPrefix
,
KSPType
kspType
,
PCType
pcType
,
MatSolverPackage
matSolverPackage
,
PetscReal
rtol
,
PetscReal
atol
,
PetscInt
maxIt
)
{
KSPSetType
(
ksp
,
kspType
);
KSPSetTolerances
(
ksp
,
rtol
,
atol
,
PETSC_DEFAULT
,
maxIt
);
...
...
@@ -267,8 +268,23 @@ namespace AMDiS {
PC
pc
;
KSPGetPC
(
ksp
,
&
pc
);
PCSetType
(
pc
,
pcType
);
if
(
matSolverPackage
!=
PETSC_NULL
)
PCFactorSetMatSolverPackage
(
pc
,
matSolverPackage
);
PCSetFromOptions
(
pc
);
}
void
setSolver
(
KSP
ksp
,
const
char
*
kspPrefix
,
KSPType
kspType
,
PCType
pcType
,
PetscReal
rtol
,
PetscReal
atol
,
PetscInt
maxIt
)
{
setSolverWithLu
(
ksp
,
kspPrefix
,
kspType
,
pcType
,
PETSC_NULL
,
rtol
,
atol
,
maxIt
);
}
}
}
AMDiS/src/parallel/PetscHelper.h
View file @
48245c1c
...
...
@@ -91,6 +91,15 @@ namespace AMDiS {
*/
void
matNestConvert
(
Mat
matNest
,
Mat
&
mat
);
void
setSolverWithLu
(
KSP
ksp
,
const
char
*
kspPrefix
,
KSPType
kspType
,
PCType
pcType
,
MatSolverPackage
matSolverPackage
,
PetscReal
rtol
=
PETSC_DEFAULT
,
PetscReal
atol
=
PETSC_DEFAULT
,
PetscInt
maxIt
=
PETSC_DEFAULT
);
void
setSolver
(
KSP
ksp
,
const
char
*
kspPrefix
,
KSPType
kspType
,
...
...
AMDiS/src/parallel/PetscSolverNavierStokes.cc
View file @
48245c1c
...
...
@@ -47,6 +47,9 @@ namespace AMDiS {
pressureComponent
(
-
1
),
pressureNullSpace
(
true
),
useOldInitialGuess
(
false
),
velocitySolutionMode
(
0
),
massSolutionMode
(
0
),
laplaceSolutionMode
(
0
),
massMatrixSolver
(
NULL
),
laplaceMatrixSolver
(
NULL
),
nu
(
NULL
),
...
...
@@ -65,6 +68,15 @@ namespace AMDiS {
Parameters
::
get
(
initFileStr
+
"->navierstokes->use old initial guess"
,
useOldInitialGuess
);
Parameters
::
get
(
initFileStr
+
"->navierstokes->velocity solver"
,
velocitySolutionMode
);
Parameters
::
get
(
initFileStr
+
"->navierstokes->mass solver"
,
massSolutionMode
);
Parameters
::
get
(
initFileStr
+
"->navierstokes->laplace solver"
,
laplaceSolutionMode
);
}
...
...
@@ -138,7 +150,19 @@ namespace AMDiS {
KSP
kspSchur
=
subKsp
[
1
];
PetscFree
(
subKsp
);
petsc_helper
::
setSolver
(
kspVelocity
,
""
,
KSPRICHARDSON
,
PCHYPRE
,
0.0
,
1e-14
,
1
);
switch
(
velocitySolutionMode
)
{
case
0
:
petsc_helper
::
setSolver
(
kspVelocity
,
""
,
KSPRICHARDSON
,
PCHYPRE
,
0.0
,
1e-14
,
1
);
break
;
case
1
:
petsc_helper
::
setSolverWithLu
(
kspVelocity
,
""
,
KSPPREONLY
,
PCLU
,
MATSOLVERMUMPS
,
0.0
,
1e-14
,
1
);
break
;
default:
ERROR_EXIT
(
"No velocity solution mode %d available!
\n
"
,
velocitySolutionMode
);
}
KSPSetType
(
kspSchur
,
KSPPREONLY
);
PC
pcSub
;
...
...
@@ -155,12 +179,17 @@ namespace AMDiS {
const
FiniteElemSpace
*
pressureFeSpace
=
componentSpaces
[
pressureComponent
];
DOFMatrix
massMatrix
(
pressureFeSpace
,
pressureFeSpace
);
Operator
massOp
(
pressureFeSpace
,
pressureFeSpace
);
if
(
!
phase
)
massOp
.
addTerm
(
new
Simple_ZOT
);
else
massOp
.
addTerm
(
new
VecAtQP_ZOT
(
phase
,
&
idFct
));
massMatrix
.
assembleOperator
(
massOp
);
{
Operator
massOp
(
pressureFeSpace
,
pressureFeSpace
);
ZeroOrderTerm
*
massTerm
=
NULL
;
if
(
!
phase
)
massTerm
=
new
Simple_ZOT
;
else
massTerm
=
new
VecAtQP_ZOT
(
phase
,
&
idFct
);
massOp
.
addTerm
(
massTerm
);
massMatrix
.
assembleOperator
(
massOp
);
delete
massTerm
;
}
massMatrixSolver
=
createSubSolver
(
pressureComponent
,
"mass_"
);
massMatrixSolver
->
fillPetscMatrix
(
&
massMatrix
);
...
...
@@ -168,12 +197,17 @@ namespace AMDiS {
// === Laplace matrix solver ===
DOFMatrix
laplaceMatrix
(
pressureFeSpace
,
pressureFeSpace
);
Operator
laplaceOp
(
pressureFeSpace
,
pressureFeSpace
);
if
(
!
phase
)
laplaceOp
.
addTerm
(
new
Simple_SOT
);
else
laplaceOp
.
addTerm
(
new
VecAtQP_SOT
(
phase
,
&
idFct
));
laplaceMatrix
.
assembleOperator
(
laplaceOp
);
{
Operator
laplaceOp
(
pressureFeSpace
,
pressureFeSpace
);
SecondOrderTerm
*
laplaceTerm
=
NULL
;
if
(
!
phase
)
laplaceTerm
=
new
Simple_SOT
;
else
laplaceTerm
=
new
VecAtQP_SOT
(
phase
,
&
idFct
);
laplaceOp
.
addTerm
(
laplaceTerm
);
laplaceMatrix
.
assembleOperator
(
laplaceOp
);
delete
laplaceTerm
;
}
laplaceMatrixSolver
=
createSubSolver
(
pressureComponent
,
string
(
"laplace_"
));
laplaceMatrixSolver
->
fillPetscMatrix
(
&
laplaceMatrix
);
...
...
@@ -190,47 +224,60 @@ namespace AMDiS {
vz
.
interpol
(
solution
->
getDOFVector
(
2
));
DOFMatrix
conDifMatrix
(
pressureFeSpace
,
pressureFeSpace
);
Operator
conDifOp
(
pressureFeSpace
,
pressureFeSpace
);
if
(
!
phase
)
{
MSG
(
"INIT WITHOUT PHASE!
\n
"
);
Simple_ZOT
*
conDif0
=
new
Simple_ZOT
(
*
invTau
);
conDifOp
.
addTerm
(
conDif0
);
Simple_SOT
*
conDif1
=
new
Simple_SOT
(
*
nu
);
conDifOp
.
addTerm
(
conDif1
);
VecAtQP_FOT
*
conDif2
=
new
VecAtQP_FOT
(
&
vx
,
&
idFct
,
0
);
conDifOp
.
addTerm
(
conDif2
,
GRD_PHI
);
VecAtQP_FOT
*
conDif3
=
new
VecAtQP_FOT
(
&
vy
,
&
idFct
,
1
);
conDifOp
.
addTerm
(
conDif3
,
GRD_PHI
);
if
(
dim
==
3
)
{
VecAtQP_FOT
*
conDif4
=
new
VecAtQP_FOT
(
&
vz
,
&
idFct
,
2
);
conDifOp
.
addTerm
(
conDif4
,
GRD_PHI
);
}
}
else
{
MSG
(
"INIT WITH PHASE!
\n
"
);
vp
.
interpol
(
phase
);
VecAtQP_ZOT
*
conDif0
=
new
VecAtQP_ZOT
(
&
vp
,
new
MultConstFct
(
*
invTau
));
conDifOp
.
addTerm
(
conDif0
);
VecAtQP_SOT
*
conDif1
=
new
VecAtQP_SOT
(
&
vp
,
new
MultConstFct
(
*
nu
));
conDifOp
.
addTerm
(
conDif1
);
Vec2AtQP_FOT
*
conDif2
=
new
Vec2AtQP_FOT
(
&
vx
,
&
vp
,
new
Multiplier3
(),
0
);
conDifOp
.
addTerm
(
conDif2
,
GRD_PHI
);
{
Operator
conDifOp
(
pressureFeSpace
,
pressureFeSpace
);
ZeroOrderTerm
*
conDif0
=
NULL
;
SecondOrderTerm
*
conDif1
=
NULL
;
FirstOrderTerm
*
conDif2
=
NULL
,
*
conDif3
=
NULL
,
*
conDif4
=
NULL
;
if
(
!
phase
)
{
MSG
(
"INIT WITHOUT PHASE!
\n
"
);
conDif0
=
new
Simple_ZOT
(
*
invTau
);
conDifOp
.
addTerm
(
conDif0
);
conDif1
=
new
Simple_SOT
(
*
nu
);
conDifOp
.
addTerm
(
conDif1
);
conDif2
=
new
VecAtQP_FOT
(
&
vx
,
&
idFct
,
0
);
conDifOp
.
addTerm
(
conDif2
,
GRD_PHI
);
conDif3
=
new
VecAtQP_FOT
(
&
vy
,
&
idFct
,
1
);
conDifOp
.
addTerm
(
conDif3
,
GRD_PHI
);
if
(
dim
==
3
)
{
conDif4
=
new
VecAtQP_FOT
(
&
vz
,
&
idFct
,
2
);
conDifOp
.
addTerm
(
conDif4
,
GRD_PHI
);
}
}
else
{
MSG
(
"INIT WITH PHASE: %e %e
\n
"
,
*
nu
,
*
invTau
);
vp
.
interpol
(
phase
);
conDif0
=
new
VecAtQP_ZOT
(
&
vp
,
new
MultConstFct
(
*
invTau
));
conDifOp
.
addTerm
(
conDif0
);
conDif1
=
new
VecAtQP_SOT
(
&
vp
,
new
MultConstFct
(
*
nu
));
conDifOp
.
addTerm
(
conDif1
);
conDif2
=
new
Vec2AtQP_FOT
(
&
vx
,
&
vp
,
new
Multiplier3
(),
0
);
conDifOp
.
addTerm
(
conDif2
,
GRD_PHI
);
conDif3
=
new
Vec2AtQP_FOT
(
&
vy
,
&
vp
,
new
Multiplier3
(),
1
);
conDifOp
.
addTerm
(
conDif3
,
GRD_PHI
);
if
(
dim
==
3
)
{
conDif4
=
new
Vec2AtQP_FOT
(
&
vz
,
&
vp
,
new
Multiplier3
(),
2
);
conDifOp
.
addTerm
(
conDif4
,
GRD_PHI
);
}
}
Vec2AtQP_FOT
*
conDif3
=
new
Vec2AtQP_FOT
(
&
vy
,
&
vp
,
new
Multiplier3
(),
1
);
conDifOp
.
addTerm
(
conDif3
,
GRD_PHI
);
conDifMatrix
.
assembleOperator
(
conDifOp
);
if
(
dim
==
3
)
{
Vec2AtQP_FOT
*
conDif4
=
new
Vec2AtQP_FOT
(
&
vz
,
&
vp
,
new
Multiplier3
(),
2
);
conDifOp
.
addTerm
(
conDif4
,
GRD_PHI
);
}
delete
conDif0
;
delete
conDif1
;
delete
conDif2
;
delete
conDif3
;
if
(
dim
==
3
)
delete
conDif4
;
}
conDifMatrix
.
assembleOperator
(
conDifOp
);
conDifMatrixSolver
=
createSubSolver
(
pressureComponent
,
"condif_"
);
conDifMatrixSolver
->
fillPetscMatrix
(
&
conDifMatrix
);
...
...
@@ -241,11 +288,52 @@ namespace AMDiS {
matShellContext
.
kspLaplace
=
laplaceMatrixSolver
->
getSolver
();
matShellContext
.
matConDif
=
conDifMatrixSolver
->
getMatInterior
();
petsc_helper
::
setSolver
(
matShellContext
.
kspMass
,
"mass_"
,
KSPCG
,
PCJACOBI
,
0.0
,
1e-14
,
2
);
petsc_helper
::
setSolver
(
matShellContext
.
kspLaplace
,
"laplace_"
,
KSPRICHARDSON
,
PCHYPRE
,
0.0
,
1e-14
,
1
);
switch
(
massSolutionMode
)
{
case
0
:
petsc_helper
::
setSolver
(
matShellContext
.
kspMass
,
"mass_"
,
KSPCG
,
PCJACOBI
,
0.0
,
1e-14
,
2
);
break
;
case
1
:
petsc_helper
::
setSolverWithLu
(
matShellContext
.
kspMass
,
"mass_"
,
KSPRICHARDSON
,
PCLU
,
MATSOLVERMUMPS
,
0.0
,
1e-14
,
1
);
break
;
default:
ERROR_EXIT
(
"No mass solution mode %d available!
\n
"
,
massSolutionMode
);
}
switch
(
laplaceSolutionMode
)
{
case
0
:
petsc_helper
::
setSolver
(
matShellContext
.
kspLaplace
,
"laplace_"
,
KSPRICHARDSON
,
PCHYPRE
,
0.0
,
1e-14
,
1
);
break
;
case
1
:
petsc_helper
::
setSolverWithLu
(
matShellContext
.
kspLaplace
,
"mass_"
,
KSPRICHARDSON
,
PCLU
,
MATSOLVERMUMPS
,
0.0
,
1e-14
,
1
);
break
;
default:
ERROR_EXIT
(
"No laplace solution mode %d available!
\n
"
,
laplaceSolutionMode
);
}
setConstantNullSpace
(
matShellContext
.
kspLaplace
);
}
void
PetscSolverNavierStokes
::
exitPreconditioner
(
PC
pc
)
{
FUNCNAME
(
"PetscSolverNavierStokes::exitPreconditioner()"
);
massMatrixSolver
->
destroyMatrixData
();
laplaceMatrixSolver
->
destroyMatrixData
();
conDifMatrixSolver
->
destroyMatrixData
();
delete
massMatrixSolver
;
massMatrixSolver
=
NULL
;
delete
laplaceMatrixSolver
;
laplaceMatrixSolver
=
NULL
;
delete
conDifMatrixSolver
;
conDifMatrixSolver
=
NULL
;
}
}
AMDiS/src/parallel/PetscSolverNavierStokes.h
View file @
48245c1c
...
...
@@ -142,6 +142,8 @@ namespace AMDiS {
void
initPreconditioner
(
PC
pc
);
void
exitPreconditioner
(
PC
pc
);
private:
int
pressureComponent
;
...
...
@@ -150,6 +152,15 @@ namespace AMDiS {
/// If true, old solution is used for initial guess in solver phase.
bool
useOldInitialGuess
;
/// 0: approximate solve 1: direct solver
int
velocitySolutionMode
;
/// 0: approximate solve 1: direct solver
int
massSolutionMode
;
/// 0: approximate solve 1: direct solver
int
laplaceSolutionMode
;
PetscSolver
*
massMatrixSolver
,
*
laplaceMatrixSolver
,
*
conDifMatrixSolver
;
NavierStokesSchurData
matShellContext
;
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment