Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
Aland, Sebastian
amdis
Commits
a6bacc8a
Commit
a6bacc8a
authored
Feb 09, 2012
by
Thomas Witkowski
Browse files
Work on repartitioning with higher order finite elements.
parent
79bcb0ed
Changes
19
Hide whitespace changes
Inline
Side-by-side
AMDiS/src/Containers.h
0 → 100644
View file @
a6bacc8a
// ============================================================================
// == ==
// == 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 Containers.h */
#ifndef AMDIS_CONTAINERS_H
#define AMDIS_CONTAINERS_H
#include <map>
#include <set>
namespace
AMDiS
{
/*
template<typename Key, typename T>
class amdis_map : public std::map<Key, T>
{
public:
inline void insert(Key &key, T &v)
{
std::map<Key, T>::insert(std::pair<Key, T>(key, v));
}
inline T& operator[](const Key& key)
{
return std::map<Key, T>::operator[](key);
}
};
*/
}
#endif
AMDiS/src/Element.cc
View file @
a6bacc8a
...
...
@@ -627,12 +627,13 @@ namespace AMDiS {
void
Element
::
getAllDofs
(
const
FiniteElemSpace
*
feSpace
,
BoundaryObject
bound
,
DofContainer
&
dofs
)
DofContainer
&
dofs
,
bool
baseDofPtr
)
{
getNodeDofs
(
feSpace
,
bound
,
dofs
);
getNodeDofs
(
feSpace
,
bound
,
dofs
,
baseDofPtr
);
if
(
feSpace
->
getBasisFcts
()
->
getDegree
()
>
1
)
getHigherOrderDofs
(
feSpace
,
bound
,
dofs
);
getHigherOrderDofs
(
feSpace
,
bound
,
dofs
,
baseDofPtr
);
}
}
AMDiS/src/Element.h
View file @
a6bacc8a
...
...
@@ -82,10 +82,8 @@ namespace AMDiS {
///
void
deleteElementDOFs
();
/** \brief
* Clone this Element and return a reference to it. Because also the DOFs
* are cloned, \ref Mesh::serializedDOfs must be used.
*/
/// Clone this Element and return a reference to it. Because also the DOFs
/// are cloned, \ref Mesh::serializedDOfs must be used.
Element
*
cloneWithDOFs
();
/** \name getting methods
...
...
@@ -111,10 +109,8 @@ namespace AMDiS {
return
child
[
i
];
}
/** \brief
* Returns true if Element is a leaf element (\ref child[0] == NULL), returns
* false otherwise.
*/
/// Returns true if Element is a leaf element (\ref child[0] == NULL), returns
/// false otherwise.
inline
const
bool
isLeaf
()
const
{
return
(
child
[
0
]
==
NULL
);
...
...
@@ -155,10 +151,8 @@ namespace AMDiS {
dofValid
=
b
;
}
/** \brief
* Returns \ref elementData's error estimation, if Element is a leaf element
* and has leaf data.
*/
/// Returns \ref elementData's error estimation, if Element is a leaf element
/// and has leaf data.
inline
double
getEstimation
(
int
row
)
const
{
if
(
isLeaf
())
{
...
...
@@ -172,10 +166,8 @@ namespace AMDiS {
return
0.0
;
}
/** \brief
* Returns Element's coarsening error estimation, if Element is a leaf
* element and if it has leaf data and if this leaf data are coarsenable.
*/
/// Returns Element's coarsening error estimation, if Element is a leaf
/// element and if it has leaf data and if this leaf data are coarsenable.
inline
double
getCoarseningEstimation
(
int
row
)
{
if
(
isLeaf
())
{
...
...
@@ -195,10 +187,8 @@ namespace AMDiS {
/// Returns local vertex number of the j-th vertex of the i-th edge
virtual
int
getVertexOfEdge
(
int
i
,
int
j
)
const
=
0
;
/** \brief
* Returns local vertex number of the vertexIndex-th vertex of the
* positionIndex-th part of type position (vertex, edge, face)
*/
/// Returns local vertex number of the vertexIndex-th vertex of the
/// positionIndex-th part of type position (vertex, edge, face)
virtual
int
getVertexOfPosition
(
GeoIndex
position
,
int
positionIndex
,
int
vertexIndex
)
const
=
0
;
...
...
@@ -215,10 +205,6 @@ namespace AMDiS {
///
virtual
DofFace
getFace
(
int
localFaceIndex
)
const
=
0
;
/// Returns either vertex, edge or face of the element.
/* template<typename T> */
/* T getObject(int index); */
/// Returns the number of parts of type i in this element
virtual
int
getGeo
(
GeoIndex
i
)
const
=
0
;
...
...
@@ -271,10 +257,8 @@ namespace AMDiS {
elementData
=
ed
;
}
/** \brief
* Sets \ref newCoord of Element. Needed by refinement, if Element has a
* boundary edge on a curved boundary.
*/
/// Sets \ref newCoord of Element. Needed by refinement, if Element has a
/// boundary edge on a curved boundary.
inline
void
setNewCoord
(
WorldVector
<
double
>*
coord
)
{
newCoord
=
coord
;
...
...
@@ -292,10 +276,8 @@ namespace AMDiS {
dof
[
pos
]
=
p
;
}
/** \brief
* Checks whether Element is a leaf element and whether it has leaf data.
* If the checks don't fail, leaf data's error estimation is set to est.
*/
/// Checks whether Element is a leaf element and whether it has leaf data.
/// If the checks don't fail, leaf data's error estimation is set to est.
inline
void
setEstimation
(
double
est
,
int
row
)
{
FUNCNAME
(
"Element::setEstimation()"
);
...
...
@@ -312,10 +294,8 @@ namespace AMDiS {
}
}
/** \brief
* Sets Element's coarsening error estimation, if Element is a leaf element
* and if it has leaf data and if this leaf data are coarsenable.
*/
/// Sets Element's coarsening error estimation, if Element is a leaf element
/// and if it has leaf data and if this leaf data are coarsenable.
inline
void
setCoarseningEstimation
(
double
est
,
int
row
)
{
if
(
isLeaf
())
{
...
...
@@ -363,18 +343,19 @@ namespace AMDiS {
virtual
const
FixVec
<
int
,
WORLD
>&
sortFaceIndices
(
int
face
,
FixVec
<
int
,
WORLD
>
*
vec
)
const
=
0
;
/// Returns a copy of itself. Needed by Mesh to create Elements by a prototype.
/// Returns a copy of itself. Needed by Mesh to create Elements by
/// a prototype.
virtual
Element
*
clone
()
=
0
;
/** \brief
* Returns which side of child[childnr] corresponds to side sidenr of this
* Element. If the child has no corresponding side, the return value is negative.
*/
/// Returns which side of child[childnr] corresponds to side sidenr of
/// this Element. If the child has no corresponding side, the return value
/// is negative.
virtual
int
getSideOfChild
(
int
childnr
,
int
sidenr
,
int
elType
=
0
)
const
=
0
;
/** \brief
* Generalization of \ref getSideOfChild to arbitrary subObject. Thus, e.g., in 3d
* we can ask for the local id of a verte, edge or face on the elements children.
* Generalization of \ref getSideOfChild to arbitrary subObject. Thus,
* e.g., in 3d we can ask for the local id of a verte, edge or face
* on the elements children.
*
* \param[in] childnr Either 0 or 1 for the left or right children.
* \param[in] subObj Defines whether we ask for VERTEX, EDGE or FACE.
...
...
@@ -384,12 +365,12 @@ namespace AMDiS {
virtual
int
getSubObjOfChild
(
int
childnr
,
GeoIndex
subObj
,
int
ithObj
,
int
elType
=
0
)
const
=
0
;
/
** \brie
f
* Returns which vertex of elements parent corresponds to the vertexnr of
* the element,
i
f
the
element is the childnr-th child of the parent
.
* If the vertex is the ner vertex at the refinement edge, -1 is returned.
*/
virtual
int
getVertexOfParent
(
int
childnr
,
int
vertexnr
,
int
elType
=
0
)
const
=
0
;
/
// Returns which vertex of elements parent corresponds to the vertexnr o
f
/// the element, if the element is the childnr-th child of the parent.
/// If the vertex
i
s
the
ner vertex at the refinement edge, -1 is returned
.
virtual
int
getVertexOfParent
(
int
childnr
,
int
vertexnr
,
int
elType
=
0
)
const
=
0
;
/// Returns whether Element is a Line
virtual
bool
isLine
()
const
=
0
;
...
...
@@ -416,14 +397,19 @@ namespace AMDiS {
* nodes alonge this vertex/edge/face are assembled and put together to
* a list.
*
* \param[in] feSpace FE space which is used to get the dofs.
* \param[in] bound Defines the vertex/edge/face of the element on
* which all vertex dofs are assembled.
* \param[out] dofs List of dofs, where the result is stored.
* \param[in] feSpace FE space which is used to get the dofs.
* \param[in] bound Defines the vertex/edge/face of the element on
* which all vertex dofs are assembled.
* \param[out] dofs List of dofs, where the result is stored.
* \param[in] baseDofPtr If true, the base DOF pointes are stored. Thus,
* dof* [\ref dof] of the element is inserted. If
* false, &(dof[.][n0]) is put to the result vector,
* with n0 beging the number of predofs.
*/
virtual
void
getNodeDofs
(
const
FiniteElemSpace
*
feSpace
,
BoundaryObject
bound
,
DofContainer
&
dofs
)
const
=
0
;
DofContainer
&
dofs
,
bool
baseDofPtr
=
false
)
const
=
0
;
/** \brief
* Traverses a vertex/edge/face of a given element (this includes also all
...
...
@@ -431,18 +417,26 @@ namespace AMDiS {
* to higher order basis functions alonge this vertex/edge/face are
* assembled and put together to a list.
*
* \param[in] feSpace FE space which is used to get the dofs.
* \param[in] bound Defines the edge/face of the element on which
* all non vertex dofs are assembled.
* \param[out] dofs All dofs are put to this dof list.
* \param[in] feSpace FE space which is used to get the dofs.
* \param[in] bound Defines the edge/face of the element on which
* all non vertex dofs are assembled.
* \param[out] dofs All dofs are put to this dof list.
* \param[in] baseDofPtr If true, the base DOF pointes are stored. Thus,
* dof* [\ref dof] of the element is inserted. If
* false, &(dof[.][n0]) is put to the result vector,
* with n0 beging the number of predofs.
*/
virtual
void
getHigherOrderDofs
(
const
FiniteElemSpace
*
feSpace
,
BoundaryObject
bound
,
DofContainer
&
dofs
)
const
=
0
;
DofContainer
&
dofs
,
bool
baseDofPtr
=
false
)
const
=
0
;
/// Combines \ref getNodeDofs and \ref getHigherOrderDofs to one function.
/// See parameter description there.
void
getAllDofs
(
const
FiniteElemSpace
*
feSpace
,
BoundaryObject
bound
,
DofContainer
&
dofs
);
DofContainer
&
dofs
,
bool
baseDofPtr
=
false
);
/** \} */
...
...
AMDiS/src/ElementDofIterator.h
View file @
a6bacc8a
...
...
@@ -82,7 +82,7 @@ namespace AMDiS {
/// Returns a pointer to the starting position of the current DOF
/// array. Makes only sence, if \ref nextStrict() is used for traverse.
inline
const
DegreeOfFreedom
*
get
Start
Dof
()
inline
const
DegreeOfFreedom
*
get
Base
Dof
()
{
return
dofs
[
node0
+
elementPos
];
}
...
...
AMDiS/src/FiniteElemSpace.cc
View file @
a6bacc8a
...
...
@@ -137,4 +137,18 @@ namespace AMDiS {
}
}
}
const
FiniteElemSpace
*
FiniteElemSpace
::
getHighest
(
vector
<
const
FiniteElemSpace
*>&
feSpaces
)
{
const
FiniteElemSpace
*
feSpace
=
feSpaces
[
0
];
for
(
unsigned
int
i
=
1
;
i
<
feSpaces
.
size
();
i
++
)
if
(
feSpaces
[
i
]
->
getBasisFcts
()
->
getDegree
()
>
feSpace
->
getBasisFcts
()
->
getDegree
())
feSpace
=
feSpaces
[
i
];
return
feSpace
;
}
}
AMDiS/src/FiniteElemSpace.h
View file @
a6bacc8a
...
...
@@ -97,6 +97,11 @@ namespace AMDiS {
int
calcMemoryUsage
();
static
void
clear
();
/// Returns for a set of FE spaces that FE space having basis functions with
/// the highest degree.
static
const
FiniteElemSpace
*
getHighest
(
vector
<
const
FiniteElemSpace
*>&
feSpaces
);
protected:
/** \brief
...
...
AMDiS/src/Line.h
View file @
a6bacc8a
...
...
@@ -173,23 +173,23 @@ namespace AMDiS {
return
"Line"
;
}
void
getNodeDofs
(
const
FiniteElemSpace
*
,
BoundaryObject
,
DofContainer
&
)
const
void
getNodeDofs
(
const
FiniteElemSpace
*
,
BoundaryObject
,
DofContainer
&
,
bool
)
const
{
FUNCNAME
(
"Line::getNodeDofs()"
);
ERROR_EXIT
(
"Not yet implemented!
\n
"
);
}
void
getHigherOrderDofs
(
const
FiniteElemSpace
*
,
BoundaryObject
,
DofContainer
&
)
const
void
getHigherOrderDofs
(
const
FiniteElemSpace
*
,
BoundaryObject
,
DofContainer
&
,
bool
)
const
{
FUNCNAME
(
"Line::getHigherOrderDofs()"
);
ERROR_EXIT
(
"Not yet implemented!
\n
"
);
}
protected:
/** \brief
* vertexOfEdge[i][j] is the local number of the j-th vertex of the i-th
* edge of this element.
*/
/// vertexOfEdge[i][j] is the local number of the j-th vertex of the i-th
/// edge of this element.
static
const
int
vertexOfEdge
[
1
][
2
];
static
const
int
sideOfChild
[
2
][
2
];
...
...
AMDiS/src/Mesh.cc
View file @
a6bacc8a
...
...
@@ -279,16 +279,10 @@ namespace AMDiS {
TEST_EXIT
(
feSpaces
.
size
()
>
0
)(
"Should not happen!
\n
"
);
// === Search for the FE space with the highest degree of polynomials. ===
// === Using this FE space ensures that deleting DOFs defined on it, ===
// === also DOFs of lower order FE spaces will be deleted correctly. ===
const
FiniteElemSpace
*
feSpace
=
feSpaces
[
0
];
for
(
unsigned
int
i
=
1
;
i
<
feSpaces
.
size
();
i
++
)
if
(
feSpaces
[
i
]
->
getBasisFcts
()
->
getDegree
()
>
feSpace
->
getBasisFcts
()
->
getDegree
())
feSpace
=
feSpaces
[
i
];
// Search for the FE space with the highest degree of polynomials. Using this
// FE space ensures that deleting DOFs defined on it, also DOFs of lower
// order FE spaces will be deleted correctly.
const
FiniteElemSpace
*
feSpace
=
FiniteElemSpace
::
getHighest
(
feSpaces
);
// === Determine to all DOFs in mesh the macro elements where the DOF ===
...
...
@@ -305,8 +299,8 @@ namespace AMDiS {
while
(
elInfo
)
{
elDofIter
.
reset
(
elInfo
->
getElement
());
do
{
dofsOwner
[
elDofIter
.
get
Start
Dof
()].
insert
(
elInfo
->
getMacroElement
());
dofsPosIndex
[
elDofIter
.
get
Start
Dof
()]
=
elDofIter
.
getPosIndex
();
dofsOwner
[
elDofIter
.
get
Base
Dof
()].
insert
(
elInfo
->
getMacroElement
());
dofsPosIndex
[
elDofIter
.
get
Base
Dof
()]
=
elDofIter
.
getPosIndex
();
}
while
(
elDofIter
.
nextStrict
());
elInfo
=
stack
.
traverseNext
(
elInfo
);
...
...
@@ -322,9 +316,9 @@ namespace AMDiS {
deque
<
MacroElement
*>
newMacroElements
;
for
(
deque
<
MacroElement
*>::
iterator
elIter
=
macroElements
.
begin
();
elIter
!=
macroElements
.
end
();
++
elIter
)
{
// If the current mesh macro element should not be deleted, i.e., it is not
a
// member of the list of macro elements to be deleted, is is inserted to
the new
// macro element list.
// If the current mesh macro element should not be deleted, i.e., it is not
//
a
member of the list of macro elements to be deleted, is is inserted to
//
the new
macro element list.
if
(
macros
.
find
(
*
elIter
)
==
macros
.
end
())
newMacroElements
.
push_back
(
*
elIter
);
}
...
...
@@ -334,15 +328,15 @@ namespace AMDiS {
macroElements
=
newMacroElements
;
// === For all macro elements to be deleted, delete them also to be
neighbours
===
// === of some other macro elements. Furtheremore, delete the
whole element
===
// === hierarchie structure of the macro element.
===
// === For all macro elements to be deleted, delete them also to be
===
// ===
neighbours
of some other macro elements. Furtheremore, delete the ===
// ===
whole element
hierarchie structure of the macro element. ===
for
(
std
::
set
<
MacroElement
*>::
iterator
macroIt
=
macros
.
begin
();
macroIt
!=
macros
.
end
();
++
macroIt
)
{
// Go through all neighbours of the macro element and remove this macro
element
// to be neighbour of some other macro element.
// Go through all neighbours of the macro element and remove this macro
//
element
to be neighbour of some other macro element.
for
(
int
i
=
0
;
i
<
getGeo
(
NEIGH
);
i
++
)
if
((
*
macroIt
)
->
getNeighbour
(
i
))
for
(
int
j
=
0
;
j
<
getGeo
(
NEIGH
);
j
++
)
...
...
@@ -361,7 +355,7 @@ namespace AMDiS {
}
// We will delete at least some element DOFs in the next but will keep the
// element object still in memory. Therefore the DOF pointer must be set
to be
// element object still in memory. Therefore the DOF pointer must be set
// invalid.
(
*
macroIt
)
->
getElement
()
->
setDofValid
(
false
);
}
...
...
AMDiS/src/MeshStructure.cc
View file @
a6bacc8a
...
...
@@ -10,15 +10,17 @@
// See also license.opensource.txt in the distribution.
#include "Debug.h"
#include "DOFVector.h"
#include "Element.h"
#include "ElementDofIterator.h"
#include "ElInfo.h"
#include "MeshStructure.h"
#include "MeshStructure_ED.h"
#include "Mesh.h"
#include "Element.h"
#include "Traverse.h"
#include "ElInfo.h"
#include "RefinementManager.h"
#include "
Debug
.h"
#include "DOFVector.h"
#include "
Traverse
.h"
namespace
AMDiS
{
...
...
@@ -416,66 +418,104 @@ namespace AMDiS {
}
void
MeshStructure
::
getMeshStructureValues
(
Mesh
*
mesh
,
int
macroElIndex
,
void
MeshStructure
::
getMeshStructureValues
(
int
macroElIndex
,
const
DOFVector
<
double
>*
vec
,
std
::
vector
<
double
>&
values
)
{
FUNCNAME
(
"MeshStructure::getMeshStructureValues()"
);
TEST_EXIT_DBG
(
mesh
)(
"No mesh defined!
\n
"
);
TEST_EXIT_DBG
(
vec
)(
"No DOFVector defined!
\n
"
);
const
FiniteElemSpace
*
feSpace
=
vec
->
getFeSpace
();
Mesh
*
mesh
=
feSpace
->
getMesh
();
bool
feSpaceHasNonVertexDofs
=
(
feSpace
->
getBasisFcts
()
->
getDegree
()
>
1
);
values
.
clear
();
ElementDofIterator
elDofIter
(
feSpace
);
TraverseStack
stack
;
ElInfo
*
elInfo
=
stack
.
traverseFirstOneMacro
(
mesh
,
macroElIndex
,
-
1
,
Mesh
::
CALL_EVERY_EL_PREORDER
);
while
(
elInfo
)
{
// For the macro element the mesh structure code stores all vertex values.
if
(
elInfo
->
getLevel
()
==
0
)
for
(
int
i
=
0
;
i
<
mesh
->
getGeo
(
VERTEX
);
i
++
)
values
.
push_back
((
*
vec
)[
elInfo
->
getElement
()
->
getDof
(
i
,
0
)]);
if
(
!
elInfo
->
getElement
()
->
isLeaf
())
values
.
push_back
((
*
vec
)[
elInfo
->
getElement
()
->
getChild
(
0
)
->
getDof
(
mesh
->
getDim
(),
0
)]);
if
(
!
elInfo
->
getElement
()
->
isLeaf
())
{
// If no leaf element store the value of the "new" DOF that is created
// by bisectioning of this element.
DegreeOfFreedom
dof0
=
elInfo
->
getElement
()
->
getChild
(
0
)
->
getDof
(
mesh
->
getDim
(),
0
);
values
.
push_back
((
*
vec
)[
dof0
]);
}
else
{
// If leaf element store all non vertex values of this element, thus
// only relevant for higher order basis functions.
if
(
feSpaceHasNonVertexDofs
)
{
elDofIter
.
reset
(
elInfo
->
getElement
());
do
{
if
(
elDofIter
.
getPosIndex
()
!=
VERTEX
)
values
.
push_back
((
*
vec
)[
elDofIter
.
getDof
()]);
}
while
(
elDofIter
.
next
());
}
}
elInfo
=
stack
.
traverseNext
(
elInfo
);
}
}
void
MeshStructure
::
setMeshStructureValues
(
Mesh
*
mesh
,
int
macroElIndex
,
void
MeshStructure
::
setMeshStructureValues
(
int
macroElIndex
,
DOFVector
<
double
>*
vec
,
const
std
::
vector
<
double
>&
values
)
{
FUNCNAME
(
"MeshStructure::setMeshStructureValues()"
);
TEST_EXIT_DBG
(
mesh
)(
"No mesh defined!
\n
"
);
TEST_EXIT_DBG
(
vec
)(
"No DOFVector defined!
\n
"
);
TEST_EXIT_DBG
(
static_cast
<
int
>
(
values
.
size
())
>=
mesh
->
getGeo
(
VERTEX
))
(
"Should not happen!
\n
"
);
const
FiniteElemSpace
*
feSpace
=
vec
->
getFeSpace
();
Mesh
*
mesh
=
feSpace
->
getMesh
();
bool
feSpaceHasNonVertexDofs
=
(
feSpace
->
getBasisFcts
()
->
getDegree
()
>
1
);
unsigned
int
counter
=
0
;
TEST_EXIT_DBG
(
static_cast
<
int
>
(
values
.
size
())
>=
mesh
->
getGeo
(
VERTEX
))
(
"Should not happen!
\n
"
);
ElementDofIterator
elDofIter
(
feSpace
);
TraverseStack
stack
;
ElInfo
*
elInfo
=
stack
.
traverseFirstOneMacro
(
mesh
,
macroElIndex
,
-
1
,
Mesh
::
CALL_EVERY_EL_PREORDER
);
while
(
elInfo
)
{
// For the macro element all vertex nodes are set first.
if
(
elInfo
->
getLevel
()
==
0
)
for
(
int
i
=
0
;
i
<
mesh
->
getGeo
(
VERTEX
);
i
++
)
(
*
vec
)[
elInfo
->
getElement
()
->
getDof
(
i
,
0
)]
=
values
[
counter
++
];
if
(
!
elInfo
->
getElement
()
->
isLeaf
())
{
// If no leaf element set the value of the "new" DOF that is created
// by bisectioning of this element.
TEST_EXIT_DBG
(
counter
<
values
.
size
())(
"Should not happen!
\n
"
);
(
*
vec
)[
elInfo
->
getElement
()
->
getChild
(
0
)
->
getDof
(
mesh
->
getDim
(),
0
)]
=
values
[
counter
++
];
}
else
{
// On leaf elements set all non vertex values (thus DOFs of higher order
// basis functions).
if
(
feSpaceHasNonVertexDofs
)
{
elDofIter
.
reset
(
elInfo
->
getElement
());
do
{
if
(
elDofIter
.
getPosIndex
()
!=
VERTEX
)
(
*
vec
)[
elDofIter
.
getDof
()]
=
values
[
counter
++
];
}
while
(
elDofIter
.
next
());
}
}
elInfo
=
stack
.
traverseNext
(
elInfo
);
}
TEST_EXIT_DBG
(
values
.
size
()
==
counter
)(
"Should not happen!
\n
"
);
}
}
AMDiS/src/MeshStructure.h
View file @
a6bacc8a
...
...
@@ -154,14 +154,32 @@ namespace AMDiS {
/// Returns true, if the given mesh structure code is equal to this one.
bool
compare
(
MeshStructure
&
other
);
void
getMeshStructureValues
(
Mesh
*
mesh
,
int
macroElIndex
,
/** \brief
* Creates a value array of a given DOFVector. This value array corresponds
* to the mesh structure code of the element and thus can easily be used
* to reconstruct the values of the DOFVector on the same element (e.g., after
* the mesh and the value array has been redistributed in parallel
* computations).
*
* \param[in] macroElIndex Index of the macro element for which the value
* structure code must be created.
* \param[in] vec DOFVector to be used for creating the value code.
* \param[out] values Resulting value structure code.
*/
void
getMeshStructureValues
(
int
macroElIndex
,
const
DOFVector
<
double
>*
vec
,
std
::
vector
<
double
>&
values
);
void
setMeshStructureValues
(
Mesh
*
mesh
,
int
macroElIndex
,
/** \brief
* Uses a value structure code, e.g. created by \ref getMeshStructureValues,
* to reconstruct the data of a DOFVector on a given macro element.
*
* \param[in] macroElIndex Macro element index the code is related to.
* \param[out] vec DOFVector that should be reconstructed.
* \param[in] values Value structure code.
*/
void
setMeshStructureValues
(
int
macroElIndex
,
DOFVector
<
double
>*
vec
,
const
std
::
vector
<
double
>&
values
);
...
...
AMDiS/src/Tetrahedron.cc
View file @
a6bacc8a
...
...
@@ -202,22 +202,27 @@ namespace AMDiS {
void
Tetrahedron
::
getNodeDofs
(
const
FiniteElemSpace
*
feSpace
,
BoundaryObject
bound
,
DofContainer
&
dofs
)
const
DofContainer
&
dofs
,
bool
baseDofPtr
)
const
{
FUNCNAME
(
"Tetrahedron::getNodeDofs()"
);
switch
(
bound
.
subObj
)
{
case
VERTEX
:
{
int
n0
=
feSpace
->
getAdmin
()
->
getNumberOfPreDofs
(
VERTEX
);
dofs
.
push_back
(
&
(
dof
[
bound
.
ithObj
][
n0
]));
if
(
baseDofPtr
)
{
dofs
.
push_back
(
dof
[
bound
.
ithObj
]);