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
2d5c10bf
Commit
2d5c10bf
authored
Feb 06, 2012
by
Thomas Witkowski
Browse files
Starting point to make FETI-DP work with mixed finite elements.
parent
988aefe0
Changes
5
Expand all
Hide whitespace changes
Inline
Side-by-side
AMDiS/src/parallel/FeSpaceMapping.h
0 → 100644
View file @
2d5c10bf
// ============================================================================
// == ==
// == AMDiS - Adaptive multidimensional simulations ==
// == ==
// == http://www.amdis-fem.org ==
// == ==
// ============================================================================
//
// Software License for AMDiS
//
// Copyright (c) 2010 Dresden University of Technology
// All rights reserved.
// Authors: Simon Vey, Thomas Witkowski et al.
//
// This file is part of AMDiS
//
// See also license.opensource.txt in the distribution.
/** \file FeSpaceMapping.h */
#include
<map>
#include
"parallel/MpiHelper.h"
#include
"parallel/ParallelTypes.h"
#ifndef AMDIS_FE_SPACE_MAPPING_H
#define AMDIS_FE_SPACE_MAPPING_H
namespace
AMDiS
{
using
namespace
std
;
class
GlobalDofMap
{
public:
GlobalDofMap
(
MPI
::
Intracomm
*
m
)
:
mpiComm
(
m
),
nRankDofs
(
0
),
nOverallDofs
(
0
),
rStartDofs
(
0
)
{}
void
clear
()
{
dofMap
.
clear
();
nRankDofs
=
0
;
nOverallDofs
=
0
;
rStartDofs
=
0
;
}
DegreeOfFreedom
operator
[](
DegreeOfFreedom
d
)
{
TEST_EXIT_DBG
(
dofMap
.
count
(
d
))(
"Should not happen!
\n
"
);
return
dofMap
[
d
];
}
void
insertRankDof
(
DegreeOfFreedom
dof0
,
DegreeOfFreedom
dof1
=
-
1
)
{
FUNCNAME
(
"GlobalDofMap::insertRankDof()"
);
TEST_EXIT_DBG
(
dofMap
.
count
(
dof0
)
==
0
)(
"Should not happen!
\n
"
);
dofMap
[
dof0
]
=
(
dof1
>=
0
?
dof1
:
nRankDofs
);
nRankDofs
++
;
}
void
insert
(
DegreeOfFreedom
dof0
,
DegreeOfFreedom
dof1
)
{
FUNCNAME
(
"GlobalDofMap::insert()"
);
TEST_EXIT_DBG
(
dofMap
.
count
(
dof0
)
==
0
)(
"Should not happen!
\n
"
);
dofMap
[
dof0
]
=
dof1
;
}
bool
isSet
(
DegreeOfFreedom
dof
)
{
return
static_cast
<
bool
>
(
dofMap
.
count
(
dof
));
}
unsigned
int
size
()
{
return
dofMap
.
size
();
}
DofMapping
&
getMap
()
{
return
dofMap
;
}
void
update
(
bool
add
=
true
)
{
nOverallDofs
=
0
;
rStartDofs
=
0
;
mpi
::
getDofNumbering
(
*
mpiComm
,
nRankDofs
,
rStartDofs
,
nOverallDofs
);
if
(
add
)
addOffset
(
rStartDofs
);
}
void
addOffset
(
int
offset
)
{
for
(
DofMapping
::
iterator
it
=
dofMap
.
begin
();
it
!=
dofMap
.
end
();
++
it
)
it
->
second
+=
offset
;
}
private:
MPI
::
Intracomm
*
mpiComm
;
///
DofMapping
dofMap
;
public:
///
int
nRankDofs
,
nOverallDofs
,
rStartDofs
;
};
template
<
typename
T
>
class
FeSpaceData
{
public:
FeSpaceData
()
{}
void
setMpiComm
(
MPI
::
Intracomm
*
m
)
{
mpiComm
=
m
;
}
T
&
operator
[](
const
FiniteElemSpace
*
feSpace
)
{
FUNCNAME
(
"FeSpaceData::operator[]()"
);
TEST_EXIT_DBG
(
data
.
count
(
feSpace
))(
"Should not happen!
\n
"
);
return
data
.
find
(
feSpace
)
->
second
;
}
void
addFeSpace
(
const
FiniteElemSpace
*
feSpace
)
{
FUNCNAME
(
"FeSpaceData::addFeSpace()"
);
if
(
data
.
count
(
feSpace
))
data
.
find
(
feSpace
)
->
second
.
clear
();
else
data
.
insert
(
make_pair
(
feSpace
,
T
(
mpiComm
)));
}
private:
MPI
::
Intracomm
*
mpiComm
;
map
<
const
FiniteElemSpace
*
,
T
>
data
;
};
}
#endif
AMDiS/src/parallel/PetscSolver.h
View file @
2d5c10bf
...
...
@@ -57,6 +57,7 @@ namespace AMDiS {
{
meshDistributor
=
m
;
mpiRank
=
meshDistributor
->
getMpiRank
();
mpiComm
=
&
(
meshDistributor
->
getMpiComm
());
}
/** \brief
...
...
@@ -85,9 +86,15 @@ namespace AMDiS {
return
0
;
}
KSP
getSolver
()
{
return
solver
;
}
KSP
getSolver
()
{
return
solver
;
}
PC
getPc
()
{
return
pc
;
}
PC
getPc
()
{
return
pc
;
}
protected:
void
printSolutionInfo
(
AdaptInfo
*
adaptInfo
,
...
...
@@ -123,6 +130,8 @@ namespace AMDiS {
int
mpiRank
;
MPI
::
Intracomm
*
mpiComm
;
/// Petsc's matrix structure.
Mat
petscMatrix
;
...
...
AMDiS/src/parallel/PetscSolverFeti.cc
View file @
2d5c10bf
This diff is collapsed.
Click to expand it.
AMDiS/src/parallel/PetscSolverFeti.h
View file @
2d5c10bf
...
...
@@ -20,7 +20,11 @@
/** \file PetscSolverFeti.h */
#include
<map>
#include
"parallel/FeSpaceMapping.h"
#include
"parallel/MpiHelper.h"
#include
"parallel/PetscSolver.h"
#include
"parallel/PetscSolverFetiStructs.h"
#include
"parallel/ParallelTypes.h"
#ifndef AMDIS_PETSC_SOLVER_FETI_H
...
...
@@ -30,100 +34,7 @@ namespace AMDiS {
using
namespace
std
;
class
PetscSolverFeti
;
/** \brief
* This structure is used when defining the MatShell operation for solving
* primal schur complement. \ref petscMultMatSchurPrimal
*/
struct
SchurPrimalData
{
/// Pointers to the matrix containing the primal variables.
Mat
*
mat_primal_primal
;
/// Coupling matrices between the primal and the B variables.
Mat
*
mat_primal_b
,
*
mat_b_primal
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vecor in the primal variables.
Vec
tmp_vec_primal
;
PetscSolverFeti
*
fetiSolver
;
};
/** \brief
* This structure is used when defining the FETI-DP operator for solving
* the system matrix reduced to the Lagrange multipliers.
* \ref petscMultMatFeti
*/
struct
FetiData
{
/// Coupling matrices between the primal and the B variables.
Mat
*
mat_primal_b
,
*
mat_b_primal
;
/// Matrix of Lagrange variables.
Mat
*
mat_lagrange
;
/// Temporal vectors on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the primal variables.
Vec
tmp_vec_primal
;
/// Temporal vector on the lagrange variables.
Vec
tmp_vec_lagrange
;
PetscSolverFeti
*
fetiSolver
;
/// Pointer to the solver of the schur complement on the primal variables.
KSP
*
ksp_schur_primal
;
};
struct
FetiDirichletPreconData
{
/// Matrix of scaled Lagrange variables.
Mat
*
mat_lagrange_scaled
;
Mat
*
mat_interior_interior
,
*
mat_duals_duals
,
*
mat_interior_duals
,
*
mat_duals_interior
;
/// Pointer to the solver for \ref PetscSolverFeti::mat_bb.
KSP
*
ksp_interior
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the dual variables.
Vec
tmp_vec_duals0
,
tmp_vec_duals1
;
/// Temporal vector on the interior variables.
Vec
tmp_vec_interior
;
};
struct
FetiLumpedPreconData
{
/// Matrix of scaled Lagrange variables.
Mat
*
mat_lagrange_scaled
;
Mat
*
mat_duals_duals
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the dual variables.
Vec
tmp_vec_duals0
,
tmp_vec_duals1
;
};
typedef
enum
{
FETI_NONE
=
0
,
FETI_DIRICHLET
=
1
,
FETI_LUMPED
=
2
}
FetiPreconditioner
;
/** \brief
* FETI-DP implementation based on PETSc.
*/
class
PetscSolverFeti
:
public
PetscSolver
...
...
@@ -238,28 +149,18 @@ namespace AMDiS {
int
nComponents
;
/// Mapping from primal DOF indices to a global index of primals.
DofMapping
globalPrimalIndex
;
/// Number of rank owned primals and global primals
int
nRankPrimals
,
nOverallPrimals
,
rStartPrimals
;
FeSpaceData
<
GlobalDofMap
>
primalDofMap
;
/// Set of DOF indices that are considered to be dual variables.
DofIndexSet
duals
;
/// Mapping from dual DOF indices to a global index of duals.
DofMapping
globalDualIndex
;
/// Stores to each dual boundary DOF the set of ranks in which the DOF
/// is contained in.
DofIndexToPartitions
boundaryDofRanks
;
FeSpaceData
<
GlobalDofMap
>
dualDofMap
;
/// Stores to each dual DOF index the index of the first Lagrange
/// constraint that is assigned to this DOF.
DofMap
ping
dofFirstLagrange
;
FeSpaceData
<
Global
DofMap
>
dofFirstLagrange
;
///
Number of rank owned Lagrange variables, number of global
///
Lagrange variables
.
int
nRankLagrange
,
nOverallLagrange
,
rStartLagrange
;
///
Stores to each dual boundary DOF the set of ranks in which the DOF
///
is contained in
.
DofIndexToPartitions
boundaryDofRanks
;
/// Index for each non primal variables to the global index of
/// B variables.
...
...
@@ -336,6 +237,7 @@ namespace AMDiS {
// Number of local nodes that are duals.
int
nLocalDuals
;
};
}
#endif
AMDiS/src/parallel/PetscSolverFetiStructs.h
0 → 100644
View file @
2d5c10bf
// ============================================================================
// == ==
// == AMDiS - Adaptive multidimensional simulations ==
// == ==
// == http://www.amdis-fem.org ==
// == ==
// ============================================================================
//
// Software License for AMDiS
//
// Copyright (c) 2010 Dresden University of Technology
// All rights reserved.
// Authors: Simon Vey, Thomas Witkowski et al.
//
// This file is part of AMDiS
//
// See also license.opensource.txt in the distribution.
/** \file PetscSolverFetiStructs.h */
#ifndef AMDIS_PETSC_SOLVER_FETI_STRUCTS_H
#define AMDIS_PETSC_SOLVER_FETI_STRUCTS_H
namespace
AMDiS
{
class
PetscSolverFeti
;
/** \brief
* This structure is used when defining the MatShell operation for solving
* primal schur complement. \ref petscMultMatSchurPrimal
*/
struct
SchurPrimalData
{
/// Pointers to the matrix containing the primal variables.
Mat
*
mat_primal_primal
;
/// Coupling matrices between the primal and the B variables.
Mat
*
mat_primal_b
,
*
mat_b_primal
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vecor in the primal variables.
Vec
tmp_vec_primal
;
PetscSolverFeti
*
fetiSolver
;
};
/** \brief
* This structure is used when defining the FETI-DP operator for solving
* the system matrix reduced to the Lagrange multipliers.
* \ref petscMultMatFeti
*/
struct
FetiData
{
/// Coupling matrices between the primal and the B variables.
Mat
*
mat_primal_b
,
*
mat_b_primal
;
/// Matrix of Lagrange variables.
Mat
*
mat_lagrange
;
/// Temporal vectors on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the primal variables.
Vec
tmp_vec_primal
;
/// Temporal vector on the lagrange variables.
Vec
tmp_vec_lagrange
;
PetscSolverFeti
*
fetiSolver
;
/// Pointer to the solver of the schur complement on the primal variables.
KSP
*
ksp_schur_primal
;
};
struct
FetiDirichletPreconData
{
/// Matrix of scaled Lagrange variables.
Mat
*
mat_lagrange_scaled
;
Mat
*
mat_interior_interior
,
*
mat_duals_duals
,
*
mat_interior_duals
,
*
mat_duals_interior
;
/// Pointer to the solver for \ref PetscSolverFeti::mat_bb.
KSP
*
ksp_interior
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the dual variables.
Vec
tmp_vec_duals0
,
tmp_vec_duals1
;
/// Temporal vector on the interior variables.
Vec
tmp_vec_interior
;
};
struct
FetiLumpedPreconData
{
/// Matrix of scaled Lagrange variables.
Mat
*
mat_lagrange_scaled
;
Mat
*
mat_duals_duals
;
/// Temporal vector on the B variables.
Vec
tmp_vec_b
;
/// Temporal vector on the dual variables.
Vec
tmp_vec_duals0
,
tmp_vec_duals1
;
};
typedef
enum
{
FETI_NONE
=
0
,
FETI_DIRICHLET
=
1
,
FETI_LUMPED
=
2
}
FetiPreconditioner
;
}
#endif
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