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
iwr
amdis
Commits
44677772
Commit
44677772
authored
Mar 23, 2009
by
Thomas Witkowski
Browse files
* small update to increase performance of assembling procedure
parent
dd8ffdf7
Changes
47
Hide whitespace changes
Inline
Side-by-side
AMDiS/compositeFEM/src/CFE_Integration.cc
View file @
44677772
...
...
@@ -279,12 +279,12 @@ namespace AMDiS {
VectorOfFixVecs
<
DimVec
<
double
>
>
&
surfVert
)
{
double
surfDet
;
int
dim
=
surfVert
[
0
].
getSize
()
-
1
;
FixVec
<
WorldVector
<
double
>
,
VERTEX
>
worldCoords
(
dim
-
1
,
NO_INIT
);
int
dim
=
surfVert
[
0
].
getSize
()
-
1
;
FixVec
<
WorldVector
<
double
>
,
VERTEX
>
worldCoords
(
dim
-
1
,
NO_INIT
);
// transform barycentric coords to world coords
for
(
int
i
=
0
;
i
<
dim
;
i
++
)
{
loc_elInfo
->
coordToWorld
(
surfVert
[
i
],
&
worldCoords
[
i
]);
for
(
int
i
=
0
;
i
<
dim
;
i
++
)
{
loc_elInfo
->
coordToWorld
(
surfVert
[
i
],
worldCoords
[
i
]);
}
// calculate determinant for surface
...
...
AMDiS/compositeFEM/src/CFE_NormAndErrorFcts.cc
View file @
44677772
...
...
@@ -17,11 +17,11 @@ namespace AMDiS {
const
double
&
fac
)
{
double
val
=
0.0
;
const
WorldVector
<
double
>
*
worldCoordsAtQP
;
WorldVector
<
double
>
worldCoordsAtQP
;
for
(
int
iq
=
0
;
iq
<
nQPts
;
++
iq
)
{
worldCoordsAtQP
=
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
NULL
);
val
+=
q
->
getWeight
(
iq
)
*
fabs
((
*
f
)(
*
worldCoordsAtQP
));
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
worldCoordsAtQP
);
val
+=
q
->
getWeight
(
iq
)
*
fabs
((
*
f
)(
worldCoordsAtQP
));
}
double
nrm
=
det
*
val
;
...
...
@@ -34,11 +34,11 @@ namespace AMDiS {
const
double
&
fac
)
{
double
val
=
0.0
;
const
WorldVector
<
double
>
*
worldCoordsAtQP
;
WorldVector
<
double
>
worldCoordsAtQP
;
for
(
int
iq
=
0
;
iq
<
nQPts
;
++
iq
)
{
worldCoordsAtQP
=
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
NULL
);
val
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
f
)(
*
worldCoordsAtQP
));
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
worldCoordsAtQP
);
val
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
f
)(
worldCoordsAtQP
));
}
double
nrm
=
det
*
val
;
...
...
@@ -51,15 +51,15 @@ namespace AMDiS {
const
double
&
fac
)
{
double
val
=
0.0
;
double
norm_grd2
;
const
WorldVector
<
double
>
*
worldCoordsAtQP
;
double
norm_grd2
=
0.0
;
WorldVector
<
double
>
worldCoordsAtQP
;
for
(
int
iq
=
0
;
iq
<
nQPts
;
++
iq
)
{
worldCoordsAtQP
=
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
NULL
);
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
worldCoordsAtQP
);
norm_grd2
=
0.0
;
for
(
int
j
=
0
;
j
<
dim
;
j
++
)
norm_grd2
+=
sqr
(((
*
grd
)(
*
worldCoordsAtQP
))[
j
]);
norm_grd2
+=
sqr
(((
*
grd
)(
worldCoordsAtQP
))[
j
]);
val
+=
q
->
getWeight
(
iq
)
*
norm_grd2
;
}
...
...
@@ -142,14 +142,14 @@ namespace AMDiS {
q
,
NULL
,
NULL
);
const
WorldVector
<
double
>
*
worldCoordsAtQP
;
WorldVector
<
double
>
worldCoordsAtQP
;
for
(
int
iq
=
0
;
iq
<
nQPts
;
++
iq
)
{
worldCoordsAtQP
=
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
NULL
);
val
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
u
)(
*
worldCoordsAtQP
)
-
uhAtQPs
[
iq
]);
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
worldCoordsAtQP
);
val
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
u
)(
worldCoordsAtQP
)
-
uhAtQPs
[
iq
]);
if
(
relErr
)
val_nrm
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
u
)(
*
worldCoordsAtQP
));
val_nrm
+=
q
->
getWeight
(
iq
)
*
sqr
((
*
u
)(
worldCoordsAtQP
));
}
double
err
=
det
*
val
;
...
...
@@ -173,23 +173,22 @@ namespace AMDiS {
q
,
NULL
,
NULL
);
const
WorldVector
<
double
>
*
worldCoordsAtQP
;
WorldVector
<
double
>
worldCoordsAtQP
;
for
(
int
iq
=
0
;
iq
<
nQPts
;
++
iq
)
{
worldCoordsAtQP
=
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
NULL
);
elInfo
->
coordToWorld
(
q
->
getLambda
(
iq
),
worldCoordsAtQP
);
norm_err_grd2
=
0.0
;
for
(
int
j
=
0
;
j
<
dim
;
++
j
)
norm_err_grd2
+=
sqr
(((
*
grdu
)(
*
worldCoordsAtQP
))[
j
]
-
grdUhAtQPs
[
iq
][
j
]);
sqr
(((
*
grdu
)(
worldCoordsAtQP
))[
j
]
-
grdUhAtQPs
[
iq
][
j
]);
val
+=
q
->
getWeight
(
iq
)
*
norm_err_grd2
;
if
(
relErr
)
{
norm_grd2
=
0.0
;
for
(
int
j
=
0
;
j
<
dim
;
++
j
)
norm_grd2
+=
sqr
(((
*
grdu
)(
*
worldCoordsAtQP
))[
j
]);
norm_grd2
+=
sqr
(((
*
grdu
)(
worldCoordsAtQP
))[
j
]);
val_nrm
+=
q
->
getWeight
(
iq
)
*
norm_grd2
;
}
...
...
AMDiS/compositeFEM/src/ElementLevelSet.cc
View file @
44677772
...
...
@@ -249,8 +249,8 @@ ElementLevelSet::calcIntersecNormal_2d(WorldVector<double> &normalVec)
// ===== Get world coordinates of intersection points. =====
WorldVector
<
double
>
sP1
;
WorldVector
<
double
>
sP2
;
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
0
],
&
sP1
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
1
],
&
sP2
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
0
],
sP1
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
1
],
sP2
);
// ===== Calculate normal vector. =====
double
norm2
=
0.0
;
...
...
@@ -326,9 +326,9 @@ ElementLevelSet::calcIntersecNormal_3d(WorldVector<double> &normalVec)
WorldVector
<
double
>
sP1
;
WorldVector
<
double
>
sP2
;
WorldVector
<
double
>
sP3
;
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
0
],
&
sP1
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
1
],
&
sP2
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
2
],
&
sP3
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
0
],
sP1
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
1
],
sP2
);
elInfo
->
coordToWorld
((
*
elIntersecPoints
)[
2
],
sP3
);
// ===== Calculate normal vector. =====
WorldVector
<
double
>
A
=
sP2
-
sP1
;
...
...
AMDiS/src/Assembler.cc
View file @
44677772
...
...
@@ -406,13 +406,11 @@ namespace AMDiS {
ZeroOrderAssembler
::
getSubAssembler
(
op
,
this
,
quad0
,
false
);
}
ElementMatrix
*
Assembler
::
initElementMatrix
(
ElementMatrix
*
elMat
,
const
ElInfo
*
rowElInfo
,
const
ElInfo
*
colElInfo
)
void
Assembler
::
initElementMatrix
(
ElementMatrix
*
elMat
,
const
ElInfo
*
rowElInfo
,
const
ElInfo
*
colElInfo
)
{
if
(
!
elMat
)
{
elMat
=
NEW
ElementMatrix
(
nRow
,
nCol
);
}
TEST_EXIT_DBG
(
elMat
)(
"No ElementMatrix!
\n
"
);
elMat
->
set
(
0.0
);
...
...
@@ -434,25 +432,17 @@ namespace AMDiS {
&
(
elMat
->
colIndices
));
}
}
return
elMat
;
}
ElementVector
*
Assembler
::
initElementVector
(
ElementVector
*
elVec
,
const
ElInfo
*
elInfo
)
void
Assembler
::
initElementVector
(
ElementVector
*
elVec
,
const
ElInfo
*
elInfo
)
{
if
(
!
elVec
)
{
elVec
=
NEW
ElementVector
(
nRow
);
}
TEST_EXIT_DBG
(
elVec
)(
"No ElementVector!
\n
"
);
elVec
->
set
(
0.0
);
DOFAdmin
*
rowAdmin
=
rowFESpace
->
getAdmin
();
Element
*
element
=
elInfo
->
getElement
();
rowFESpace
->
getBasisFcts
()
->
getLocalIndicesVec
(
element
,
rowAdmin
,
&
(
elVec
->
dofIndices
));
return
elVec
;
rowFESpace
->
getBasisFcts
()
->
getLocalIndicesVec
(
elInfo
->
getElement
(),
rowFESpace
->
getAdmin
(),
&
(
elVec
->
dofIndices
));
}
void
Assembler
::
checkQuadratures
()
...
...
AMDiS/src/Assembler.h
View file @
44677772
...
...
@@ -69,13 +69,13 @@ namespace AMDiS {
virtual
~
Assembler
()
{}
ElementMatrix
*
initElementMatrix
(
ElementMatrix
*
elMat
,
const
ElInfo
*
rowElInfo
,
const
ElInfo
*
colElInfo
=
NULL
);
void
initElementMatrix
(
ElementMatrix
*
elMat
,
const
ElInfo
*
rowElInfo
,
const
ElInfo
*
colElInfo
=
NULL
);
ElementVector
*
initElementVector
(
ElementVector
*
elVec
,
const
ElInfo
*
elInfo
);
void
initElementVector
(
ElementVector
*
elVec
,
const
ElInfo
*
elInfo
);
/** \brief
* Assembles the element matrix for the given ElInfo
...
...
@@ -230,9 +230,7 @@ namespace AMDiS {
const
ElInfo
*
smallElInfo
,
const
ElInfo
*
largeElInfo
,
ElementVector
*
vec
);
/** \brief
* Checks whether this is a new travese.
*/
/// Checks whether this is a new travese.
inline
void
checkForNewTraverse
()
{
if
(
lastTraverseId
!=
ElInfo
::
traverseId
[
omp_get_thread_num
()])
{
lastVecEl
=
lastMatEl
=
NULL
;
...
...
AMDiS/src/BoundaryManager.cc
View file @
44677772
...
...
@@ -12,6 +12,7 @@ namespace AMDiS {
BoundaryManager
::
BoundaryManager
(
const
FiniteElemSpace
*
feSpace
)
{
localBounds
.
resize
(
omp_get_overall_max_threads
());
dofIndices
.
resize
(
omp_get_overall_max_threads
());
allocatedMemoryLocalBounds
=
feSpace
->
getBasisFcts
()
->
getNumber
();
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
localBounds
.
size
());
i
++
)
{
localBounds
[
i
]
=
GET_MEMORY
(
BoundaryType
,
allocatedMemoryLocalBounds
);
...
...
@@ -23,6 +24,7 @@ namespace AMDiS {
localBCs
=
bm
.
localBCs
;
allocatedMemoryLocalBounds
=
bm
.
allocatedMemoryLocalBounds
;
localBounds
.
resize
(
bm
.
localBounds
.
size
());
dofIndices
.
resize
(
bm
.
localBounds
.
size
());
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
localBounds
.
size
());
i
++
)
{
localBounds
[
i
]
=
GET_MEMORY
(
BoundaryType
,
allocatedMemoryLocalBounds
);
}
...
...
@@ -51,29 +53,26 @@ namespace AMDiS {
void
BoundaryManager
::
fillBoundaryConditions
(
ElInfo
*
elInfo
,
DOFVectorBase
<
double
>
*
vec
)
{
// ===== fill local conditions ==============================================
const
FiniteElemSpace
*
feSpace
=
vec
->
getFESpace
();
Vector
<
DegreeOfFreedom
>
dofIndices
;
const
BasisFunction
*
basisFcts
=
feSpace
->
getBasisFcts
();
int
nBasFcts
=
basisFcts
->
getNumber
();
DOFAdmin
*
admin
=
feSpace
->
getAdmin
();
std
::
map
<
BoundaryType
,
BoundaryCondition
*>::
iterator
it
;
if
(
localBCs
.
size
()
>
0
)
{
const
FiniteElemSpace
*
feSpace
=
vec
->
getFESpace
();
Vector
<
DegreeOfFreedom
>
&
dofVec
=
dofIndices
[
omp_get_thread_num
()];
const
BasisFunction
*
basisFcts
=
feSpace
->
getBasisFcts
();
int
nBasFcts
=
basisFcts
->
getNumber
();
std
::
map
<
BoundaryType
,
BoundaryCondition
*>::
iterator
it
;
// get boundaries of all DOFs
BoundaryType
*
localBound
=
localBounds
[
omp_get_thread_num
()];
basisFcts
->
getBound
(
elInfo
,
localBound
);
// get dof indices
basisFcts
->
getLocalIndicesVec
(
elInfo
->
getElement
(),
admin
,
&
dofIndices
);
basisFcts
->
getLocalIndicesVec
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
&
dofVec
);
// apply non dirichlet boundary conditions
for
(
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
{
if
((
*
it
).
second
)
{
if
(
!
(
*
it
).
second
->
isDirichlet
())
{
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dof
Indices
[
0
],
localBound
,
nBasFcts
);
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dof
Vec
[
0
],
localBound
,
nBasFcts
);
}
}
}
...
...
@@ -82,7 +81,7 @@ namespace AMDiS {
for
(
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
{
if
((
*
it
).
second
)
{
if
((
*
it
).
second
->
isDirichlet
())
{
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dof
Indices
[
0
],
localBound
,
nBasFcts
);
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dof
Vec
[
0
],
localBound
,
nBasFcts
);
}
}
}
...
...
@@ -92,28 +91,26 @@ namespace AMDiS {
void
BoundaryManager
::
fillBoundaryConditions
(
ElInfo
*
elInfo
,
DOFMatrix
*
mat
)
{
// ===== fill local conditions ==============================================
const
FiniteElemSpace
*
feSpace
=
mat
->
getRowFESpace
();
Vector
<
DegreeOfFreedom
>
dofIndices
;
const
BasisFunction
*
basisFcts
=
feSpace
->
getBasisFcts
();
int
nBasFcts
=
basisFcts
->
getNumber
();
DOFAdmin
*
admin
=
feSpace
->
getAdmin
();
std
::
map
<
BoundaryType
,
BoundaryCondition
*>::
iterator
it
;
if
(
localBCs
.
size
()
>
0
)
{
const
FiniteElemSpace
*
feSpace
=
mat
->
getRowFESpace
();
Vector
<
DegreeOfFreedom
>
&
dofVec
=
dofIndices
[
omp_get_thread_num
()];
const
BasisFunction
*
basisFcts
=
feSpace
->
getBasisFcts
();
int
nBasFcts
=
basisFcts
->
getNumber
();
std
::
map
<
BoundaryType
,
BoundaryCondition
*>::
iterator
it
;
// get boundaries of all DOFs
BoundaryType
*
localBound
=
localBounds
[
omp_get_thread_num
()];
basisFcts
->
getBound
(
elInfo
,
localBound
);
// get dof indices
basisFcts
->
getLocalIndicesVec
(
elInfo
->
getElement
(),
admin
,
&
dofIndices
);
basisFcts
->
getLocalIndicesVec
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
&
dofVec
);
// apply non dirichlet boundary conditions
for
(
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
{
if
((
*
it
).
second
)
{
if
(
!
(
*
it
).
second
->
isDirichlet
())
{
(
*
it
).
second
->
fillBoundaryCondition
(
mat
,
elInfo
,
&
dof
Indices
[
0
],
localBound
,
nBasFcts
);
(
*
it
).
second
->
fillBoundaryCondition
(
mat
,
elInfo
,
&
dof
Vec
[
0
],
localBound
,
nBasFcts
);
}
}
}
...
...
@@ -122,7 +119,7 @@ namespace AMDiS {
for
(
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
{
if
((
*
it
).
second
)
{
if
((
*
it
).
second
->
isDirichlet
())
{
(
*
it
).
second
->
fillBoundaryCondition
(
mat
,
elInfo
,
&
dof
Indices
[
0
],
localBound
,
nBasFcts
);
(
*
it
).
second
->
fillBoundaryCondition
(
mat
,
elInfo
,
&
dof
Vec
[
0
],
localBound
,
nBasFcts
);
}
}
}
...
...
AMDiS/src/BoundaryManager.h
View file @
44677772
...
...
@@ -31,6 +31,7 @@ namespace AMDiS {
class
DOFMatrix
;
class
FiniteElemSpace
;
template
<
typename
T
>
class
Vector
;
template
<
typename
T
>
class
DOFVectorBase
;
// ============================================================================
...
...
@@ -107,16 +108,15 @@ namespace AMDiS {
};
protected:
/** \brief
* Map of managed local boundary conditions.
*/
/// Map of managed local boundary conditions.
std
::
map
<
BoundaryType
,
BoundaryCondition
*>
localBCs
;
/** \brief
* Temporary thread-safe variable for functions fillBoundaryconditions.
*/
/// Temporary thread-safe variable for functions fillBoundaryconditions.
std
::
vector
<
BoundaryType
*>
localBounds
;
/// Temporary thread-safe variable for functions fillBoundaryconditions.
std
::
vector
<
Vector
<
DegreeOfFreedom
>
>
dofIndices
;
/** \brief
* Stores the number of byte that were allocated in the constructor for
* each localBounds value. Is used to free the memory in the destructor.
...
...
AMDiS/src/DOFMatrix.cc
View file @
44677772
...
...
@@ -249,10 +249,6 @@ namespace AMDiS {
int
nRow
=
elMat
.
rowIndices
.
getSize
();
int
nCol
=
elMat
.
colIndices
.
getSize
();
// elMat.rowIndices.print();
// elMat.colIndices.print();
// elMat.print();
for
(
int
i
=
0
;
i
<
nRow
;
i
++
)
{
// for all rows of element matrix
row
=
elMat
.
rowIndices
[
i
];
BoundaryCondition
*
condition
=
...
...
@@ -433,37 +429,40 @@ namespace AMDiS {
}
void
DOFMatrix
::
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
,
Operator
*
op
)
void
DOFMatrix
::
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
)
{
FUNCNAME
(
"DOFMatrix::assemble()"
);
if
(
!
op
&&
operators
.
size
()
==
0
)
{
return
;
}
if
(
operators
.
size
()
==
0
)
return
;
Operator
*
operat
=
op
?
op
:
operators
[
0
];
operat
->
getAssembler
(
omp_get_thread_num
())
->
initElementMatrix
(
elementMatrix
,
elInfo
);
operators
[
0
]
->
getAssembler
(
omp_get_thread_num
())
->
initElementMatrix
(
elementMatrix
,
elInfo
);
if
(
op
)
{
op
->
getElementMatrix
(
elInfo
,
elementMatrix
);
}
else
{
std
::
vector
<
Operator
*>::
iterator
it
;
std
::
vector
<
double
*>::
iterator
factorIt
;
for
(
it
=
operators
.
begin
(),
factorIt
=
operatorFactor
.
begin
();
it
!=
operators
.
end
();
++
it
,
++
factorIt
)
{
std
::
vector
<
Operator
*>::
iterator
it
=
operators
.
begin
();
std
::
vector
<
double
*>::
iterator
factorIt
=
operatorFactor
.
begin
();
for
(;
it
!=
operators
.
end
();
++
it
,
++
factorIt
)
{
if
((
*
it
)
->
getNeedDualTraverse
()
==
false
)
{
(
*
it
)
->
getElementMatrix
(
elInfo
,
elementMatrix
,
*
factorIt
?
**
factorIt
:
1.0
);
(
*
it
)
->
getElementMatrix
(
elInfo
,
elementMatrix
,
*
factorIt
?
**
factorIt
:
1.0
);
}
}
}
}
addElementMatrix
(
factor
,
*
elementMatrix
,
bound
);
}
void
DOFMatrix
::
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
,
Operator
*
op
)
{
FUNCNAME
(
"DOFMatrix::assemble()"
);
TEST_EXIT_DBG
(
op
)(
"No operator!
\n
"
);
op
->
getAssembler
(
omp_get_thread_num
())
->
initElementMatrix
(
elementMatrix
,
elInfo
);
op
->
getElementMatrix
(
elInfo
,
elementMatrix
);
addElementMatrix
(
factor
,
*
elementMatrix
,
bound
);
}
void
DOFMatrix
::
assemble
(
double
factor
,
ElInfo
*
rowElInfo
,
ElInfo
*
colElInfo
,
ElInfo
*
smallElInfo
,
ElInfo
*
largeElInfo
,
...
...
AMDiS/src/DOFMatrix.h
View file @
44677772
...
...
@@ -357,9 +357,10 @@ namespace AMDiS {
* matrix should be cleared by calling clear dof matrix().
*/
void
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
,
Operator
*
op
=
NULL
);
void
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
);
void
assemble
(
double
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
,
Operator
*
op
);
void
assemble
(
double
factor
,
ElInfo
*
rowElInfo
,
ElInfo
*
colElInfo
,
...
...
AMDiS/src/DOFVector.cc
View file @
44677772
...
...
@@ -511,7 +511,6 @@ namespace AMDiS {
elInfo
=
stack
.
traverseNext
(
elInfo
);
}
}
else
{
DimVec
<
double
>
*
coords1
=
NULL
;
DimVec
<
double
>
coords2
(
feSpace
->
getMesh
()
->
getDim
(),
NO_INIT
);
DualTraverse
dualStack
;
ElInfo
*
elInfo1
,
*
elInfo2
;
...
...
@@ -534,8 +533,7 @@ namespace AMDiS {
for
(
int
i
=
0
;
i
<
nBasisFcts
;
i
++
)
{
if
(
vec
[
myLocalIndices
[
i
]]
==
0.0
)
{
coords1
=
basisFcts
->
getCoords
(
i
);
elInfo1
->
coordToWorld
(
*
coords1
,
&
worldVec
);
elInfo1
->
coordToWorld
(
*
(
basisFcts
->
getCoords
(
i
)),
worldVec
);
elInfo2
->
worldToCoord
(
worldVec
,
&
coords2
);
bool
isPositive
=
true
;
...
...
@@ -821,13 +819,14 @@ namespace AMDiS {
{
FUNCNAME
(
"DOFVector::assemble()"
);
TEST_EXIT_DBG
(
this
->
elementVector
)(
"No ElementVector!
\n
"
);
if
(
!
(
op
||
this
->
operators
.
size
()))
return
;
Operator
*
operat
=
op
?
op
:
this
->
operators
[
0
];
this
->
elementVector
=
operat
->
getAssembler
(
omp_get_thread_num
())
->
initElementVector
(
this
->
elementVector
,
elInfo
);
operat
->
getAssembler
(
omp_get_thread_num
())
->
initElementVector
(
this
->
elementVector
,
elInfo
);
if
(
op
)
{
op
->
getElementVector
(
elInfo
,
this
->
elementVector
);
...
...
@@ -862,9 +861,8 @@ namespace AMDiS {
Operator
*
operat
=
op
?
op
:
this
->
operators
[
0
];
this
->
elementVector
=
operat
->
getAssembler
(
omp_get_thread_num
())
->
initElementVector
(
this
->
elementVector
,
mainElInfo
);
operat
->
getAssembler
(
omp_get_thread_num
())
->
initElementVector
(
this
->
elementVector
,
mainElInfo
);
if
(
op
)
{
ERROR_EXIT
(
"TODO"
);
...
...
AMDiS/src/DOFVector.h
View file @
44677772
...
...
@@ -560,7 +560,7 @@ namespace AMDiS {
/** \brief
* Assignment operator between two vectors
*/
DOFVector
<
T
>&
operator
=
(
const
DOFVector
<
T
>&
);
DOFVector
<
T
>&
operator
=
(
const
DOFVector
<
T
>&
);
/** \brief
* vec[i] = v.vec[i]
...
...
AMDiS/src/DOFVector.hh
View file @
44677772
...
...
@@ -26,7 +26,6 @@ namespace AMDiS {
DOFVectorBase
<
T
>::
DOFVectorBase
(
const
FiniteElemSpace
*
f
,
std
::
string
n
)
:
feSpace
(
f
),
name
(
n
),
elementVector
(
NULL
),
boundaryManager
(
NULL
)
{
nBasFcts
=
feSpace
->
getBasisFcts
()
->
getNumber
();
...
...
@@ -45,6 +44,8 @@ namespace AMDiS {
grdTmp
[
i
]
=
NEW
DimVec
<
double
>
(
dim
,
DEFAULT_VALUE
,
0.0
);
D2Phis
[
i
]
=
NEW
DimMat
<
double
>
(
dim
,
NO_INIT
);
}
elementVector
=
NEW
ElementVector
(
this
->
nBasFcts
);
}
template
<
typename
T
>
...
...
@@ -743,8 +744,9 @@ namespace AMDiS {
DOFVector
<
T
>&
DOFVector
<
T
>::
operator
=
(
const
DOFVector
<
T
>&
rhs
)
{
feSpace
=
rhs
.
feSpace
;
this
->
nBasFcts
=
rhs
.
nBasFcts
;
vec
=
rhs
.
vec
;
this
->
elementVector
=
NULL
;
this
->
elementVector
=
new
ElementVector
(
this
->
nBasFcts
)
;
interFct
=
rhs
.
interFct
;
coarsenOperation
=
rhs
.
coarsenOperation
;
this
->
operators
=
rhs
.
operators
;
...
...
@@ -756,7 +758,7 @@ namespace AMDiS {
// boundaryManager->setDOFVector(this);
}
else
{
this
->
boundaryManager
=
NULL
;
this
->
boundaryManager
=
NULL
;
}
return
*
this
;
...
...
AMDiS/src/DataCollector.cc
View file @
44677772
...
...
@@ -220,7 +220,7 @@ namespace AMDiS {
const
DegreeOfFreedom
**
dof
=
elInfo
->
getElement
()
->
getDOF
();
DegreeOfFreedom
vertexDOF
;
WorldVector
<
double
>
vertexCoords
;
// create ElementInfo
ElementInfo
elementInfo
(
dim_
);
...
...
@@ -326,7 +326,7 @@ namespace AMDiS {
for
(
int
i
=
mesh_
->
getGeo
(
VERTEX
);
i
<
nBasisFcts_
;
i
++
)
{
DegreeOfFreedom
dofi
=
localDOFs_
[
i
];
elInfo
->
coordToWorld
(
*
basisFcts_
->
getCoords
(
i
),
vertexCoords
);
elInfo
->
coordToWorld
(
*
basisFcts_
->
getCoords
(
i
),
*
vertexCoords
);
if
((
*
interpPointInd_
)[
dofi
]
==
-
1
)
{
// mark as interpolation point
...
...
AMDiS/src/DataCollector.h
View file @
44677772
...
...
@@ -213,44 +213,28 @@ namespace AMDiS {
*/
int
nPreDofs_
;
/** \brief
* Number of vertices.
*/
/// Number of vertices.
int
nVertices_
;
/** \brief
* Number of elements.
*/
/// Number of elements.
int
nElements_
;
/** \brief
* Total number of interpolation points.
*/
/// Total number of interpolation points.
int
nInterpPoints_
;
/** \brief
* Number of connections in periodic problems.
*/
/// Number of connections in periodic problems.
int
nConnections_
;
/** \brief
* Dimension of \ref mesh
*/
/// Dimension of \ref mesh
int
dim_
;
/** \brief
* Maps internal element indices to global output indices.
*/
/// Maps internal element indices to global output indices.
std
::
map
<
int
,
int
>
outputIndices_
;
/** \brief
* Global interpolation point indexing
*/
/// Global interpolation point indexing
DOFVector
<
int
>
*
interpPointInd_
;