Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-gfe
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Sander, Oliver
dune-gfe
Commits
39d2b367
Commit
39d2b367
authored
13 years ago
by
Youett, Jonathan
Committed by
akbib@FU-BERLIN.DE
13 years ago
Browse files
Options
Downloads
Patches
Plain Diff
make class inherit from BasisGridFunction
[[Imported from SVN: r7972]]
parent
de6578f0
No related branches found
No related tags found
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dune/gfe/globalgfetestfunction.hh
+16
-20
16 additions, 20 deletions
dune/gfe/globalgfetestfunction.hh
with
16 additions
and
20 deletions
dune/gfe/globalgfetestfunction.hh
+
16
−
20
View file @
39d2b367
...
...
@@ -10,17 +10,21 @@
#include
<dune/gfe/localgfetestfunction.hh>
#include
<dune/fufem/functions/
virtual
gridfunction.hh>
#include
<dune/fufem/functions/
basis
gridfunction.hh>
/** \brief Global geodesic finite element test function.
*
* \tparam B
asis
- The global basis type.
* \tparam B - The global basis type.
* \tparam TargetSpace - The manifold that this functions takes its values in.
* \tparam CoefficientType - The coefficient vector type.
*/
template
<
class
Basis
,
class
TargetSpace
,
class
CoefficientType
>
GlobalGFETestFunction
:
public
VirtualGridFunction
<
typename
Basis
::
GridView
::
Grid
,
typename
TargetSpace
::
EmbeddedTangentVector
>
{
template
<
class
B
,
class
TargetSpace
,
class
CoefficientType
>
GlobalGFETestFunction
:
public
BasisGridFunction
<
Basis
,
CoefficientType
>
{
public:
typedef
B
Basis
;
private:
typedef
typename
Basis
::
LocalFiniteElement
LocalFiniteElement
;
typedef
typename
Basis
::
GridView
GridView
;
typedef
typename
GridView
::
template
Codim
<
0
>
::
Entity
Element
;
...
...
@@ -38,11 +42,9 @@ GlobalGFETestFunction : public VirtualGridFunction<typename Basis::GridView::Gri
enum
{
tangentDim
=
TargetSpace
::
TangentVector
::
dimension
};
public
:
//! Create global function by a global gfe test basis and the corresponding coefficient vectors
GlobalGFETestFunction
(
const
Basis
&
basis
,
const
CoefficientType
&
coefficients
)
:
basis_
(
basis
),
coefficients_
(
coefficients
)
BasisGridFunction
<
Basis
,
CoefficientType
>
(
basis
,
coefficients
)
{}
/** \brief Evaluate the function at local coordinates. */
...
...
@@ -52,33 +54,27 @@ public:
void
evaluateDerivativeLocal
(
const
Element
&
element
,
const
Dune
::
FieldVector
<
gridDim
,
ctype
>&
local
,
Dune
::
FieldMatrix
<
ctype
,
embeddedDim
,
gridDim
>&
out
)
const
;
private
:
//! The global basis
const
Basis
&
basis_
;
//! The coefficient vector for the global test function
const
CoefficientType
&
coefficients_
;
};
template
<
class
Basis
,
class
TargetSpace
,
class
CoefficientType
>
void
GlobalGFETestFunction
<
Basis
,
TargetSpace
>::
evaluateLocal
(
const
Element
&
element
,
const
Dune
::
FieldVector
<
gridDim
,
ctype
>&
local
,
EmbeddedTangentVector
&
out
)
const
{
int
numOfBasisFct
=
basis_
.
getLocalFiniteElement
(
element
).
size
();
int
numOfBasisFct
=
this
->
basis_
.
getLocalFiniteElement
(
element
).
size
();
// values of the test basis functions
std
::
vector
<
Dune
::
array
<
EmbeddedTangentVector
,
tangentDim
>
>
values
;
// create local gfe test function
basis_
.
getLocalFiniteElement
(
element
).
localBasis
().
evaluateFunction
(
local
,
values
);
this
->
basis_
.
getLocalFiniteElement
(
element
).
localBasis
().
evaluateFunction
(
local
,
values
);
// multiply values with the corresponding test coefficients and sum them up
out
=
0
;
for
(
int
i
=
0
;
i
<
values
.
size
();
i
++
)
{
int
index
=
basis_
.
index
(
element
,
i
);
int
index
=
this
->
basis_
.
index
(
element
,
i
);
for
(
int
j
=
0
;
j
<
tangentDim
;
j
++
)
{
values
[
i
][
j
]
*=
coefficients_
[
index
][
j
];
values
[
i
][
j
]
*=
this
->
coefficients_
[
index
][
j
];
out
+=
values
[
i
][
j
];
}
}
...
...
@@ -92,15 +88,15 @@ void GlobalGFETestFunction<Basis,TargetSpace>::evaluateDerivativeLocal(const Ele
std
::
vector
<
Dune
::
array
<
Dune
::
FieldMatrix
<
ctype
,
embeddedDim
,
gridDim
>
,
tangentDim
>
>
jacobians
;
// evaluate local gfe test function basis
basis_
.
getLocalFiniteElement
(
element
).
localBasis
().
evaluateJacobian
(
local
,
jacobians
);
this
->
basis_
.
getLocalFiniteElement
(
element
).
localBasis
().
evaluateJacobian
(
local
,
jacobians
);
// multiply values with the corresponding test coefficients and sum them up
out
=
0
;
for
(
int
i
=
0
;
i
<
jacobians
.
size
();
i
++
)
{
int
index
=
basis_
.
index
(
element
,
i
);
int
index
=
this
->
basis_
.
index
(
element
,
i
);
for
(
int
j
=
0
;
j
<
tangentDim
;
j
++
)
{
jacobians
[
i
][
j
]
*=
coefficients_
[
index
][
j
];
jacobians
[
i
][
j
]
*=
this
->
coefficients_
[
index
][
j
];
out
+=
jacobians
[
i
][
j
];
}
}
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
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!
Save comment
Cancel
Please
register
or
sign in
to comment