Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-fufem
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Requirements
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Locked files
Build
Pipelines
Jobs
Pipeline schedules
Test cases
Artifacts
Deploy
Releases
Container registry
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Code review analytics
Issue analytics
Insights
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
Ansgar Burchardt
dune-fufem
Commits
88e1d449
Commit
88e1d449
authored
8 years ago
by
Carsten Gräser
Browse files
Options
Downloads
Patches
Plain Diff
Add functional assembler for Dune::Functions bases
parent
adc06619
Branches
Branches containing commit
Tags
Tags containing commit
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
dune/fufem/assemblers/CMakeLists.txt
+1
-0
1 addition, 0 deletions
dune/fufem/assemblers/CMakeLists.txt
dune/fufem/assemblers/dunefunctionsfunctionalassembler.hh
+89
-0
89 additions, 0 deletions
dune/fufem/assemblers/dunefunctionsfunctionalassembler.hh
with
90 additions
and
0 deletions
dune/fufem/assemblers/CMakeLists.txt
+
1
−
0
View file @
88e1d449
...
...
@@ -5,6 +5,7 @@ install(FILES
boundaryfunctionalassembler.hh
boundaryoperatorassembler.hh
dunefunctionsoperatorassembler.hh
dunefunctionsfunctionalassembler.hh
functionalassembler.hh
integraloperatorassembler.hh
intersectionoperatorassembler.hh
...
...
This diff is collapsed.
Click to expand it.
dune/fufem/assemblers/dunefunctionsfunctionalassembler.hh
0 → 100644
+
89
−
0
View file @
88e1d449
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_FUFEM_ASSEMBERS_DUNEFUNCTIONSFUNCTIONALASSEMBLER_HH
#define DUNE_FUFEM_ASSEMBERS_DUNEFUNCTIONSFUNCTIONALASSEMBLER_HH
#include
<dune/istl/matrix.hh>
#include
<dune/istl/matrixindexset.hh>
#include
<dune/matrix-vector/axpy.hh>
#include
"dune/fufem/functionspacebases/functionspacebasis.hh"
namespace
Dune
{
namespace
Fufem
{
//! Generic global assembler for functionals on a gridview
template
<
class
TrialBasis
>
class
DuneFunctionsFunctionalAssembler
{
using
GridView
=
typename
TrialBasis
::
GridView
;
public:
//! create assembler for grid
DuneFunctionsFunctionalAssembler
(
const
TrialBasis
&
tBasis
)
:
trialBasis_
(
tBasis
)
{}
template
<
class
VectorBackend
,
class
LocalAssembler
>
void
assembleBulkEntries
(
VectorBackend
&&
vectorBackend
,
LocalAssembler
&&
localAssembler
)
const
{
auto
trialLocalView
=
trialBasis_
.
localView
();
auto
trialLocalIndexSet
=
trialBasis_
.
localIndexSet
();
using
Field
=
std
::
decay_t
<
decltype
(
vectorBackend
(
trialLocalIndexSet
.
index
(
0
)))
>
;
using
LocalVector
=
Dune
::
BlockVector
<
Dune
::
FieldVector
<
Field
,
1
>>
;
auto
localVector
=
LocalVector
(
trialLocalView
.
maxSize
());
for
(
const
auto
&
element
:
elements
(
trialBasis_
.
gridView
()))
{
trialLocalView
.
bind
(
element
);
trialLocalIndexSet
.
bind
(
trialLocalView
);
localAssembler
(
element
,
localVector
,
trialLocalView
);
// Add element stiffness matrix onto the global stiffness matrix
for
(
size_t
localRow
=
0
;
localRow
<
trialLocalIndexSet
.
size
();
++
localRow
)
{
auto
row
=
trialLocalIndexSet
.
index
(
localRow
);
vectorBackend
(
row
)
+=
localVector
[
localRow
];
}
}
}
template
<
class
VectorBackend
,
class
LocalAssembler
>
void
assembleBulk
(
VectorBackend
&&
vectorBackend
,
LocalAssembler
&&
localAssembler
)
const
{
vectorBackend
.
resize
(
trialBasis_
);
vectorBackend
.
vector
()
=
0.0
;
assembleBulkEntries
(
vectorBackend
,
std
::
forward
<
LocalAssembler
>
(
localAssembler
));
}
protected
:
const
TrialBasis
&
trialBasis_
;
};
/**
* \brief Create DuneFunctionsFunctionalAssembler
*/
template
<
class
TrialBasis
>
auto
duneFunctionsFunctionalAssembler
(
const
TrialBasis
&
trialBasis
)
{
return
DuneFunctionsFunctionalAssembler
<
TrialBasis
>
(
trialBasis
);
}
}
// namespace Fufem
}
// namespace Dune
#endif // DUNE_FUFEM_ASSEMBERS_DUNEFUNCTIONSFUNCTIONALASSEMBLER_HH
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