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
9bbe8339
Commit
9bbe8339
authored
Jan 25, 2012
by
Thomas Witkowski
Browse files
evalAtPoint
parent
7856d906
Changes
9
Show whitespace changes
Inline
Side-by-side
AMDiS/src/BoundaryManager.cc
View file @
9bbe8339
...
@@ -71,16 +71,20 @@ namespace AMDiS {
...
@@ -71,16 +71,20 @@ namespace AMDiS {
basisFcts
->
getBound
(
elInfo
,
localBound
);
basisFcts
->
getBound
(
elInfo
,
localBound
);
// get dof indices
// get dof indices
basisFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
dofVec
);
basisFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
dofVec
);
// apply non dirichlet boundary conditions
// apply non dirichlet boundary conditions
for
(
BoundaryIndexMap
::
iterator
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
for
(
BoundaryIndexMap
::
iterator
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
if
((
*
it
).
second
&&
!
(
*
it
).
second
->
isDirichlet
())
if
((
*
it
).
second
&&
!
(
*
it
).
second
->
isDirichlet
())
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dofVec
[
0
],
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dofVec
[
0
],
localBound
,
nBasFcts
);
localBound
,
nBasFcts
);
// apply dirichlet boundary conditions
// apply dirichlet boundary conditions
for
(
BoundaryIndexMap
::
iterator
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
for
(
BoundaryIndexMap
::
iterator
it
=
localBCs
.
begin
();
it
!=
localBCs
.
end
();
++
it
)
if
((
*
it
).
second
&&
(
*
it
).
second
->
isDirichlet
())
if
((
*
it
).
second
&&
(
*
it
).
second
->
isDirichlet
())
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dofVec
[
0
],
(
*
it
).
second
->
fillBoundaryCondition
(
vec
,
elInfo
,
&
dofVec
[
0
],
localBound
,
nBasFcts
);
localBound
,
nBasFcts
);
...
@@ -101,7 +105,7 @@ namespace AMDiS {
...
@@ -101,7 +105,7 @@ namespace AMDiS {
// get boundaries of all DOFs
// get boundaries of all DOFs
basisFcts
->
getBound
(
elInfo
,
localBound
);
basisFcts
->
getBound
(
elInfo
,
localBound
);
// get
dof
indices
// get
DOF
indices
basisFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
dofVec
);
basisFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
dofVec
);
// apply non dirichlet boundary conditions
// apply non dirichlet boundary conditions
...
...
AMDiS/src/BoundaryManager.h
View file @
9bbe8339
...
@@ -137,12 +137,12 @@ namespace AMDiS {
...
@@ -137,12 +137,12 @@ namespace AMDiS {
/** \brief
/** \brief
* For every boundary id we store here all possible boundary object (although
* For every boundary id we store here all possible boundary object (although
* it's not clear if it is meaningful to have different boundary conditions on
the
* it's not clear if it is meaningful to have different boundary conditions on
* same boundary id).
*
the
same boundary id).
*
*
* We have to use this global variable, because the mesh traverse interface
does
* We have to use this global variable, because the mesh traverse interface
* not provide more information about traversed boundaries at elements
than the
*
does
not provide more information about traversed boundaries at elements
* boundary id.
*
than the
boundary id.
*
*
* TODO: Change interface such that mesh traverse returns the boundary objects
* TODO: Change interface such that mesh traverse returns the boundary objects
* directly and we can remove this global variable. The biggest problem will be
* directly and we can remove this global variable. The biggest problem will be
...
...
AMDiS/src/DOFVector.cc
View file @
9bbe8339
...
@@ -68,9 +68,10 @@ namespace AMDiS {
...
@@ -68,9 +68,10 @@ namespace AMDiS {
template
<
>
template
<
>
const
double
&
DOFVector
<
double
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
,
double
*
values
)
const
const
double
DOFVector
<
double
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
)
const
{
{
FUNCNAME
(
"DOFVector<double>::evalAt
Coords
()"
);
FUNCNAME
(
"DOFVector<double>::evalAt
Point
()"
);
Mesh
*
mesh
=
feSpace
->
getMesh
();
Mesh
*
mesh
=
feSpace
->
getMesh
();
const
BasisFunction
*
basFcts
=
feSpace
->
getBasisFcts
();
const
BasisFunction
*
basFcts
=
feSpace
->
getBasisFcts
();
...
@@ -78,7 +79,7 @@ namespace AMDiS {
...
@@ -78,7 +79,7 @@ namespace AMDiS {
int
dim
=
mesh
->
getDim
();
int
dim
=
mesh
->
getDim
();
int
nBasFcts
=
basFcts
->
getNumber
();
int
nBasFcts
=
basFcts
->
getNumber
();
DegreeOfFreedom
*
localIndices
=
new
DegreeOfFreedom
[
nBasFcts
]
;
std
::
vector
<
DegreeOfFreedom
>
localIndices
(
nBasFcts
)
;
DimVec
<
double
>
lambda
(
dim
,
NO_INIT
);
DimVec
<
double
>
lambda
(
dim
,
NO_INIT
);
ElInfo
*
elInfo
=
mesh
->
createNewElInfo
();
ElInfo
*
elInfo
=
mesh
->
createNewElInfo
();
...
@@ -86,7 +87,8 @@ namespace AMDiS {
...
@@ -86,7 +87,8 @@ namespace AMDiS {
bool
inside
=
false
;
bool
inside
=
false
;
if
(
oldElInfo
&&
oldElInfo
->
getMacroElement
())
{
if
(
oldElInfo
&&
oldElInfo
->
getMacroElement
())
{
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
oldElInfo
->
getMacroElement
(),
NULL
,
NULL
);
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
oldElInfo
->
getMacroElement
(),
NULL
,
NULL
);
delete
oldElInfo
;
delete
oldElInfo
;
}
else
}
else
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
NULL
,
NULL
,
NULL
);
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
NULL
,
NULL
,
NULL
);
...
@@ -95,7 +97,8 @@ namespace AMDiS {
...
@@ -95,7 +97,8 @@ namespace AMDiS {
oldElInfo
=
elInfo
;
oldElInfo
=
elInfo
;
if
(
inside
)
{
if
(
inside
)
{
basFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
localIndices
);
basFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
localIndices
);
ElementVector
uh
(
lambda
.
size
());
ElementVector
uh
(
lambda
.
size
());
for
(
int
i
=
0
;
i
<
lambda
.
size
();
i
++
)
for
(
int
i
=
0
;
i
<
lambda
.
size
();
i
++
)
uh
[
i
]
=
operator
[](
localIndices
[
i
]);
uh
[
i
]
=
operator
[](
localIndices
[
i
]);
...
@@ -103,20 +106,19 @@ namespace AMDiS {
...
@@ -103,20 +106,19 @@ namespace AMDiS {
}
else
}
else
throw
(
std
::
runtime_error
(
"Can not eval DOFVector at point p, because point is outside geometry."
));
throw
(
std
::
runtime_error
(
"Can not eval DOFVector at point p, because point is outside geometry."
));
delete
[]
localIndices
;
if
(
oldElInfo
==
NULL
)
if
(
oldElInfo
==
NULL
)
delete
elInfo
;
delete
elInfo
;
if
(
values
!=
NULL
)
*
values
=
value
;
return
value
;
return
value
;
}
;
}
template
<
>
template
<
>
const
WorldVector
<
double
>&
DOFVector
<
WorldVector
<
double
>
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
,
WorldVector
<
double
>*
values
)
const
const
WorldVector
<
double
>
DOFVector
<
WorldVector
<
double
>
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
)
const
{
{
FUNCNAME
(
"DOFVector<double>::evalAt
Coords
()"
);
FUNCNAME
(
"DOFVector<double>::evalAt
Point
()"
);
Mesh
*
mesh
=
feSpace
->
getMesh
();
Mesh
*
mesh
=
feSpace
->
getMesh
();
const
BasisFunction
*
basFcts
=
feSpace
->
getBasisFcts
();
const
BasisFunction
*
basFcts
=
feSpace
->
getBasisFcts
();
...
@@ -124,18 +126,17 @@ namespace AMDiS {
...
@@ -124,18 +126,17 @@ namespace AMDiS {
int
dim
=
mesh
->
getDim
();
int
dim
=
mesh
->
getDim
();
int
nBasFcts
=
basFcts
->
getNumber
();
int
nBasFcts
=
basFcts
->
getNumber
();
DegreeOfFreedom
*
localIndices
=
new
DegreeOfFreedom
[
nBasFcts
]
;
vector
<
DegreeOfFreedom
>
localIndices
(
nBasFcts
)
;
DimVec
<
double
>
lambda
(
dim
,
NO_INIT
);
DimVec
<
double
>
lambda
(
dim
,
NO_INIT
);
ElInfo
*
elInfo
=
mesh
->
createNewElInfo
();
ElInfo
*
elInfo
=
mesh
->
createNewElInfo
();
static
WorldVector
<
double
>
Values
(
DEFAULT_VALUE
,
0.0
);
WorldVector
<
double
>
value
(
DEFAULT_VALUE
,
0.0
);
WorldVector
<
double
>
*
val
=
(
NULL
!=
values
)
?
values
:
&
Values
;
bool
inside
=
false
;
bool
inside
=
false
;
if
(
oldElInfo
&&
oldElInfo
->
getMacroElement
())
{
if
(
oldElInfo
&&
oldElInfo
->
getMacroElement
())
{
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
oldElInfo
->
getMacroElement
(),
NULL
,
NULL
);
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
oldElInfo
->
getMacroElement
(),
NULL
,
NULL
);
delete
oldElInfo
;
delete
oldElInfo
;
}
else
}
else
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
NULL
,
NULL
,
NULL
);
inside
=
mesh
->
findElInfoAtPoint
(
p
,
elInfo
,
lambda
,
NULL
,
NULL
,
NULL
);
...
@@ -144,20 +145,20 @@ namespace AMDiS {
...
@@ -144,20 +145,20 @@ namespace AMDiS {
oldElInfo
=
elInfo
;
oldElInfo
=
elInfo
;
if
(
inside
)
{
if
(
inside
)
{
basFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
localIndices
);
basFcts
->
getLocalIndices
(
elInfo
->
getElement
(),
feSpace
->
getAdmin
(),
localIndices
);
mtl
::
dense_vector
<
WorldVector
<
double
>
>
uh
(
lambda
.
size
());
mtl
::
dense_vector
<
WorldVector
<
double
>
>
uh
(
lambda
.
size
());
for
(
int
i
=
0
;
i
<
lambda
.
size
();
i
++
)
for
(
int
i
=
0
;
i
<
lambda
.
size
();
i
++
)
uh
[
i
]
=
operator
[](
localIndices
[
i
]);
uh
[
i
]
=
operator
[](
localIndices
[
i
]);
*
val
=
basFcts
->
evalUh
(
lambda
,
uh
);
val
ue
=
basFcts
->
evalUh
(
lambda
,
uh
);
}
else
}
else
throw
(
std
::
runtime_error
(
"Can not eval DOFVector at point p, because point is outside geometry."
));
throw
(
std
::
runtime_error
(
"Can not eval DOFVector at point p, because point is outside geometry."
));
delete
[]
localIndices
;
if
(
oldElInfo
==
NULL
)
if
(
oldElInfo
==
NULL
)
delete
elInfo
;
delete
elInfo
;
return
((
*
val
))
;
return
val
ue
;
}
;
}
template
<
>
template
<
>
...
...
AMDiS/src/DOFVector.h
View file @
9bbe8339
...
@@ -65,17 +65,13 @@ namespace AMDiS {
...
@@ -65,17 +65,13 @@ namespace AMDiS {
virtual
~
DOFVectorBase
();
virtual
~
DOFVectorBase
();
/** \brief
/// For the given element, this function returns an array of all DOFs of this
* For the given element, this function returns an array of all DOFs of this
/// DOFVector that are defined on this element.
* DOFVector that are defined on this element.
*/
virtual
void
getLocalVector
(
const
Element
*
el
,
virtual
void
getLocalVector
(
const
Element
*
el
,
mtl
::
dense_vector
<
T
>&
localVec
)
const
;
mtl
::
dense_vector
<
T
>&
localVec
)
const
;
/** \brief
/// Evaluates the DOF vector at a set of quadrature points defined on the
* Evaluates the DOF vector at a set of quadrature points defined on the
/// given element.
* given element.
*/
void
getVecAtQPs
(
const
ElInfo
*
elInfo
,
void
getVecAtQPs
(
const
ElInfo
*
elInfo
,
const
Quadrature
*
quad
,
const
Quadrature
*
quad
,
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
...
@@ -87,10 +83,8 @@ namespace AMDiS {
...
@@ -87,10 +83,8 @@ namespace AMDiS {
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
mtl
::
dense_vector
<
T
>&
vecAtQPs
)
const
;
mtl
::
dense_vector
<
T
>&
vecAtQPs
)
const
;
/** \brief
/// Evaluates the gradient of a DOF vector at a set of quadrature points
* Evaluates the gradient of a DOF vector at a set of quadrature points defined on the
/// defined on the given element.
* given element.
*/
void
getGrdAtQPs
(
const
ElInfo
*
elInfo
,
void
getGrdAtQPs
(
const
ElInfo
*
elInfo
,
const
Quadrature
*
quad
,
const
Quadrature
*
quad
,
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
...
@@ -102,10 +96,8 @@ namespace AMDiS {
...
@@ -102,10 +96,8 @@ namespace AMDiS {
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
mtl
::
dense_vector
<
typename
GradientType
<
T
>::
type
>
&
grdAtQPs
)
const
;
mtl
::
dense_vector
<
typename
GradientType
<
T
>::
type
>
&
grdAtQPs
)
const
;
/** \brief
/// Evaluates the comp'th component of the derivative of a DOF vector at a
* Evaluates the comp'th component of the derivative of a DOF vector at a set of quadrature points defined on the
/// set of quadrature points defined on the given element.
* given element.
*/
void
getDerivativeAtQPs
(
const
ElInfo
*
elInfo
,
void
getDerivativeAtQPs
(
const
ElInfo
*
elInfo
,
const
Quadrature
*
quad
,
const
Quadrature
*
quad
,
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
...
@@ -119,10 +111,8 @@ namespace AMDiS {
...
@@ -119,10 +111,8 @@ namespace AMDiS {
int
comp
,
int
comp
,
mtl
::
dense_vector
<
T
>
&
derivativeAtQPs
)
const
;
mtl
::
dense_vector
<
T
>
&
derivativeAtQPs
)
const
;
/** \brief
/// Evaluates the jacobian of a DOF vector at a set of quadrature points
* Evaluates the jacobian of a DOF vector at a set of quadrature points defined on the
/// defined on the given element.
* given element.
*/
void
getD2AtQPs
(
const
ElInfo
*
elInfo
,
void
getD2AtQPs
(
const
ElInfo
*
elInfo
,
const
Quadrature
*
quad
,
const
Quadrature
*
quad
,
const
FastQuadrature
*
quadFast
,
const
FastQuadrature
*
quadFast
,
...
@@ -134,10 +124,8 @@ namespace AMDiS {
...
@@ -134,10 +124,8 @@ namespace AMDiS {
return
feSpace
;
return
feSpace
;
}
}
/** \brief
/// Assembles the element vector for the given ellement and adds the element
* Assembles the element vector for the given ellement and adds the
/// matrix to the current DOF vector.
* element matrix to the current DOF vector.
*/
void
assemble
(
T
factor
,
ElInfo
*
elInfo
,
void
assemble
(
T
factor
,
ElInfo
*
elInfo
,
const
BoundaryType
*
bound
,
const
BoundaryType
*
bound
,
Operator
*
op
=
NULL
);
Operator
*
op
=
NULL
);
...
@@ -154,11 +142,9 @@ namespace AMDiS {
...
@@ -154,11 +142,9 @@ namespace AMDiS {
ElInfo
*
elInfo
,
ElInfo
*
elInfo
,
bool
add
=
true
);
bool
add
=
true
);
/* \brief
/// That function must be called after the matrix assembling has been
* That function must be called after the matrix assembling has been finished.
/// finished. This makes it possible to start some cleanup or matrix data
* This makes it possible to start some cleanup or matrix data compressing
/// compressing procedures.
* procedures.
*/
void
finishAssembling
();
void
finishAssembling
();
inline
void
addOperator
(
Operator
*
op
,
inline
void
addOperator
(
Operator
*
op
,
...
@@ -391,10 +377,8 @@ namespace AMDiS {
...
@@ -391,10 +377,8 @@ namespace AMDiS {
return
vec
.
end
();
return
vec
.
end
();
}
}
/** \brief
/// Used by DOFAdmin to compress this DOFVector. Implementation of
* Used by DOFAdmin to compress this DOFVector. Implementation of
/// DOFIndexedBase::compress()
* DOFIndexedBase::compress()
*/
virtual
void
compressDOFIndexed
(
int
first
,
int
last
,
virtual
void
compressDOFIndexed
(
int
first
,
int
last
,
std
::
vector
<
DegreeOfFreedom
>
&
newDof
);
std
::
vector
<
DegreeOfFreedom
>
&
newDof
);
...
@@ -473,10 +457,8 @@ namespace AMDiS {
...
@@ -473,10 +457,8 @@ namespace AMDiS {
/// Calculates Integral of this DOFVector
/// Calculates Integral of this DOFVector
double
Int
(
Quadrature
*
q
=
NULL
)
const
;
double
Int
(
Quadrature
*
q
=
NULL
)
const
;
/** \brief
/// Calculates Integral of this DOFVector over parts of the domain
* Calculates Integral of this DOFVector over parts of the domain
/// boundary, indicated by boundaryType. Implemented for DOFVector<double>
* boundary, indicated by boundaryType. Implemented for DOFVector<double>
**/
T
IntOnBoundary
(
BoundaryType
boundary
,
Quadrature
*
q
=
NULL
)
const
T
IntOnBoundary
(
BoundaryType
boundary
,
Quadrature
*
q
=
NULL
)
const
{
{
FUNCNAME
(
"DOFVector::IntOnBoundary())"
);
FUNCNAME
(
"DOFVector::IntOnBoundary())"
);
...
@@ -484,11 +466,9 @@ namespace AMDiS {
...
@@ -484,11 +466,9 @@ namespace AMDiS {
return
0.0
;
return
0.0
;
}
}
/** \brief
/// Calculates Integral of this DOFVector times normal vector over parts
* Calculates Integral of this DOFVector times normal vector over parts
/// of the domain boundary, indicated by boundaryType. Implemented for
* of the domain boundary, indicated by boundaryType. Implemented for
/// DOFVector<WorldVector<double> >
* DOFVector<WorldVector<double> >
**/
double
IntOnBoundaryNormal
(
BoundaryType
boundary
,
Quadrature
*
q
=
NULL
)
const
double
IntOnBoundaryNormal
(
BoundaryType
boundary
,
Quadrature
*
q
=
NULL
)
const
{
{
FUNCNAME
(
"DOFVector::IntOnBoundaryNormal())"
);
FUNCNAME
(
"DOFVector::IntOnBoundaryNormal())"
);
...
@@ -581,25 +561,21 @@ namespace AMDiS {
...
@@ -581,25 +561,21 @@ namespace AMDiS {
///
///
int
calcMemoryUsage
()
const
;
int
calcMemoryUsage
()
const
;
/** \brief
/// Computes the coefficients of the interpolant of the function fct and
* Computes the coefficients of the interpolant of the function fct and
/// stores these in the DOFVector
* stores these in the DOFVector
*/
void
interpol
(
AbstractFunction
<
T
,
WorldVector
<
double
>
>
*
fct
);
void
interpol
(
AbstractFunction
<
T
,
WorldVector
<
double
>
>
*
fct
);
void
interpol
(
DOFVector
<
T
>
*
v
,
double
factor
=
1.0
);
void
interpol
(
DOFVector
<
T
>
*
v
,
double
factor
=
1.0
);
/** eval DOFVector at given point p. If oldElInfo != NULL the search for the element, where p is inside,
/// Eval DOFVector at given point p. If oldElInfo != NULL the search for
* starts from oldElInfo.
/// the element, where p is inside, starts from oldElInfo. implemented for:
* implemented for: double, WorldVector< double >
/// double, WorldVector<double>
*/
inline
const
T
evalAtPoint
(
WorldVector
<
double
>
&
p
,
inline
const
T
&
evalAtPoint
(
ElInfo
*
oldElInfo
=
NULL
)
const
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
=
NULL
,
T
*
value
=
NULL
)
const
{
{
FUNCNAME
(
"DOFVector::evalAtPoint())"
);
FUNCNAME
(
"DOFVector::evalAtPoint())"
);
TEST_EXIT
(
false
)(
"Please implement your evaluation
\n
"
);
TEST_EXIT
(
false
)(
"Please implement your evaluation
\n
"
);
return
*
value
;
}
}
/// Writes the data of the DOFVector to an output stream.
/// Writes the data of the DOFVector to an output stream.
...
@@ -619,13 +595,11 @@ namespace AMDiS {
...
@@ -619,13 +595,11 @@ namespace AMDiS {
in
.
read
(
reinterpret_cast
<
char
*>
(
&
(
vec
[
0
])),
size
*
sizeof
(
T
));
in
.
read
(
reinterpret_cast
<
char
*>
(
&
(
vec
[
0
])),
size
*
sizeof
(
T
));
}
}
// DOFVector<WorldVector<T> > *getGradient(DOFVector<WorldVector<T> >*) const;
DOFVector
<
typename
GradientType
<
T
>::
type
>*
DOFVector
<
typename
GradientType
<
T
>::
type
>*
getGradient
(
DOFVector
<
typename
GradientType
<
T
>::
type
>
*
grad
)
const
;
getGradient
(
DOFVector
<
typename
GradientType
<
T
>::
type
>
*
grad
)
const
;
WorldVector
<
DOFVector
<
T
>*>
*
getGradient
(
WorldVector
<
DOFVector
<
T
>*>
*
grad
)
const
;
WorldVector
<
DOFVector
<
T
>*>
*
getGradient
(
WorldVector
<
DOFVector
<
T
>*>
*
grad
)
const
;
// DOFVector<WorldVector<T> >* getRecoveryGradient(DOFVector<WorldVector<T> >*) const;
DOFVector
<
typename
GradientType
<
T
>::
type
>*
DOFVector
<
typename
GradientType
<
T
>::
type
>*
getRecoveryGradient
(
DOFVector
<
typename
GradientType
<
T
>::
type
>
*
grad
)
const
;
getRecoveryGradient
(
DOFVector
<
typename
GradientType
<
T
>::
type
>
*
grad
)
const
;
...
@@ -654,12 +628,13 @@ namespace AMDiS {
...
@@ -654,12 +628,13 @@ namespace AMDiS {
BoundaryType
boundaryType
,
Quadrature
*
q
)
const
;
BoundaryType
boundaryType
,
Quadrature
*
q
)
const
;
template
<
>
template
<
>
const
double
&
DOFVector
<
double
>::
evalAtPoint
(
const
double
DOFVector
<
double
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
,
double
*
value
)
const
;
ElInfo
*
oldElInfo
)
const
;
template
<
>
template
<
>
const
WorldVector
<
double
>&
DOFVector
<
WorldVector
<
double
>
>::
evalAtPoint
(
const
WorldVector
<
double
>
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
,
WorldVector
<
double
>*
value
)
const
;
DOFVector
<
WorldVector
<
double
>
>::
evalAtPoint
(
WorldVector
<
double
>
&
p
,
ElInfo
*
oldElInfo
)
const
;
template
<
>
template
<
>
void
DOFVector
<
double
>::
refineInterpol
(
RCNeighbourList
&
,
int
);
void
DOFVector
<
double
>::
refineInterpol
(
RCNeighbourList
&
,
int
);
...
@@ -686,10 +661,8 @@ namespace AMDiS {
...
@@ -686,10 +661,8 @@ namespace AMDiS {
public
DOFContainer
public
DOFContainer
{
{
public:
public:
/** \brief
/// Calls constructor of DOFVector<DegreeOfFreedom> and registers itself
* Calls constructor of DOFVector<DegreeOfFreedom> and registers itself
/// as DOFContainer at DOFAdmin
* as DOFContainer at DOFAdmin
*/
DOFVectorDOF
(
const
FiniteElemSpace
*
feSpace_
,
std
::
string
name_
)
DOFVectorDOF
(
const
FiniteElemSpace
*
feSpace_
,
std
::
string
name_
)
:
DOFVector
<
DegreeOfFreedom
>
(
feSpace_
,
name_
)
:
DOFVector
<
DegreeOfFreedom
>
(
feSpace_
,
name_
)
{
{
...
@@ -702,10 +675,8 @@ namespace AMDiS {
...
@@ -702,10 +675,8 @@ namespace AMDiS {
feSpace
->
getAdmin
()
->
removeDOFContainer
(
this
);
feSpace
->
getAdmin
()
->
removeDOFContainer
(
this
);
}
}
/** \brief
/// Implements DOFContainer::operator[]() by calling
* Implements DOFContainer::operator[]() by calling
/// DOFVector<DegreeOfFreedom>::operator[]()
* DOFVector<DegreeOfFreedom>::operator[]()
*/
DegreeOfFreedom
&
operator
[](
DegreeOfFreedom
i
)
DegreeOfFreedom
&
operator
[](
DegreeOfFreedom
i
)
{
{
return
DOFVector
<
DegreeOfFreedom
>::
operator
[](
i
);
return
DOFVector
<
DegreeOfFreedom
>::
operator
[](
i
);
...
@@ -869,20 +840,22 @@ namespace AMDiS {
...
@@ -869,20 +840,22 @@ namespace AMDiS {
std
::
vector
<
DOFVector
<
double
>*>
*
res
);
std
::
vector
<
DOFVector
<
double
>*>
*
res
);
/** \brief
/** \brief
* Computes the integral of a function that includes two different DOFVectors.
This
* Computes the integral of a function that includes two different DOFVectors.
* function works also for the case that the DOFVectors are defined on
two different
*
This
function works also for the case that the DOFVectors are defined on
* meshes.
*
two different
meshes.
*/
*/
double
integrate
(
const
DOFVector
<
double
>
&
vec1
,
double
integrate
(
const
DOFVector
<
double
>
&
vec1
,
const
DOFVector
<
double
>
&
vec2
,
const
DOFVector
<
double
>
&
vec2
,
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
/// Computes the integral of a function with two DOFVectors defined on the same mesh.
/// Computes the integral of a function with two DOFVectors defined on the
/// same mesh.
double
intSingleMesh
(
const
DOFVector
<
double
>
&
vec1
,
double
intSingleMesh
(
const
DOFVector
<
double
>
&
vec1
,
const
DOFVector
<
double
>
&
vec2
,
const
DOFVector
<
double
>
&
vec2
,
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
/// Computes the integral of a function with two DOFVectors defined on different meshes.
/// Computes the integral of a function with two DOFVectors defined on
/// different meshes.
double
intDualMesh
(
const
DOFVector
<
double
>
&
vec1
,
double
intDualMesh
(
const
DOFVector
<
double
>
&
vec1
,
const
DOFVector
<
double
>
&
vec2
,
const
DOFVector
<
double
>
&
vec2
,
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
BinaryAbstractFunction
<
double
,
double
,
double
>
*
fct
);
...
...
AMDiS/src/Mesh.cc
View file @
9bbe8339
...
@@ -589,7 +589,7 @@ namespace AMDiS {
...
@@ -589,7 +589,7 @@ namespace AMDiS {
{
{
const
DOFAdmin
*
localAdmin
=
NULL
;
const
DOFAdmin
*
localAdmin
=
NULL
;
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
admin
.
size
()
)
;
i
++
)
{
for
(
unsigned
int
i
=
0
;
i
<
admin
.
size
();
i
++
)
{
if
(
admin
[
i
]
->
getNumberOfDofs
(
VERTEX
))
{
if
(
admin
[
i
]
->
getNumberOfDofs
(
VERTEX
))
{
if
(
!
localAdmin
)
if
(
!
localAdmin
)
localAdmin
=
admin
[
i
];
localAdmin
=
admin
[
i
];
...
...
AMDiS/src/Mesh.h
View file @
9bbe8339
...
@@ -223,7 +223,8 @@ namespace AMDiS {
...
@@ -223,7 +223,8 @@ namespace AMDiS {
/// Returns a DOFAdmin which at least manages vertex DOFs
/// Returns a DOFAdmin which at least manages vertex DOFs
const
DOFAdmin
*
getVertexAdmin
()
const
;
const
DOFAdmin
*
getVertexAdmin
()
const
;
/// Allocates a array of DOF pointers. The array holds one pointer for each node.
/// Allocates an array of DOF pointers. The array holds one pointer for
/// each node.
DegreeOfFreedom
**
createDofPtrs
();
DegreeOfFreedom
**
createDofPtrs
();
/// Returns \ref preserveCoarseDOFs of the mesh
/// Returns \ref preserveCoarseDOFs of the mesh
...
@@ -850,7 +851,8 @@ namespace AMDiS {
...
@@ -850,7 +851,8 @@ namespace AMDiS {
*/
*/
Element
*
elementPrototype
;
Element
*
elementPrototype
;
/// Prototype for leaf data. Used for creation of new leaf data while refinement.
/// Prototype for leaf data. Used for creation of new leaf data while
/// refinement.
ElementData
*
elementDataPrototype
;
ElementData
*
elementDataPrototype
;
/// Used for enumeration of all mesh elements
/// Used for enumeration of all mesh elements
...
@@ -862,23 +864,20 @@ namespace AMDiS {
...
@@ -862,23 +864,20 @@ namespace AMDiS {
/// Map of managed periodic vertex associations.
/// Map of managed periodic vertex associations.
map
<
BoundaryType
,
VertexVector
*>
periodicAssociations
;
map
<
BoundaryType
,
VertexVector
*>
periodicAssociations
;