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
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package registry
Model registry
Operate
Environments
Terraform modules
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
fufem
dune-fufem
Commits
9d4742e7
Commit
9d4742e7
authored
2 months ago
by
Carsten Gräser
Browse files
Options
Downloads
Patches
Plain Diff
[constraints] Add utitily function for computing Dirichlet constraints
parent
288f78e2
Branches
Branches containing commit
Tags
Tags containing commit
1 merge request
!258
Add constraints support for dune-functions basis
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
dune/fufem/constraints/CMakeLists.txt
+1
-0
1 addition, 0 deletions
dune/fufem/constraints/CMakeLists.txt
dune/fufem/constraints/boundaryconstraints.hh
+92
-0
92 additions, 0 deletions
dune/fufem/constraints/boundaryconstraints.hh
with
93 additions
and
0 deletions
dune/fufem/constraints/CMakeLists.txt
+
1
−
0
View file @
9d4742e7
install
(
FILES
affineconstraints.hh
boundaryconstraints.hh
continuityconstraints.hh
DESTINATION
${
CMAKE_INSTALL_INCLUDEDIR
}
/dune/fufem/constraints
)
This diff is collapsed.
Click to expand it.
dune/fufem/constraints/boundaryconstraints.hh
0 → 100644
+
92
−
0
View file @
9d4742e7
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_FUFEM_CONSTRAINTS_BOUNDARYCONSTRAINTS_HH
#define DUNE_FUFEM_CONSTRAINTS_BOUNDARYCONSTRAINTS_HH
#include
<variant>
#include
<dune/functions/functionspacebases/interpolate.hh>
#include
<dune/fufem/backends/traversal.hh>
#include
<dune/fufem/boundarypatch.hh>
#include
<dune/fufem/functiontools/boundarydofs.hh>
#include
<dune/fufem/constraints/affineconstraints.hh>
namespace
Dune
::
Fufem
{
/**
* \brief Compute constraints for essential boundary conditions
* \ingroup Constraints
*
* \tparam Basis Type of the global basis
* \tparam Function Type of the function to interpolate
*
* \param constraints An initialized AffineConstraints object
* \param basis A global basis
* \param f The Dirichlet values given as a function
* \param boundaryPatch The BoundaryPath where the boundary condition should be applied
*
* This will interpolate the function `f` on the `boundaryPatch`
* add constraints for all corresponding DOFs to `constraints`.
*
* This method assumes that the `constraints` has already been initialized.
* If `constraints` already contains constraint DOFs, then interpolation
* and the newly added constraints are reduced to the previously unconstrained
* DOFs.
*
* \warning This function is still experimental and the interface may change.
*/
template
<
class
BV
,
class
V
,
class
MI
,
class
C
,
class
Basis
,
class
Function
>
auto
computeBoundaryConstraints
(
AffineConstraints
<
BV
,
V
,
MI
,
C
>&
constraints
,
const
Basis
&
basis
,
Function
&&
f
,
const
BoundaryPatch
<
typename
Basis
::
GridView
>&
boundaryPatch
)
{
auto
&&
isConstrained
=
Dune
::
Functions
::
istlVectorBackend
(
constraints
.
isConstrained
());
auto
isNewConstrained
=
constraints
.
isConstrained
();
Dune
::
Fufem
::
markBoundaryPatchDofs
(
boundaryPatch
,
basis
,
isNewConstrained
);
Dune
::
Fufem
::
recursiveForEachVectorEntry
(
isNewConstrained
,
[
&
](
auto
&&
isNewConstrained_i
,
const
auto
&
i
)
{
isNewConstrained_i
=
isNewConstrained_i
and
not
isConstrained
[
i
];
isConstrained
[
i
]
=
isNewConstrained_i
or
isConstrained
[
i
];
});
Dune
::
Functions
::
interpolate
(
basis
,
constraints
.
constantTerm
(),
f
,
isNewConstrained
);
constraints
.
resolveDependencies
();
}
/**
* \brief Compute constraints for essential boundary conditions
* \ingroup Constraints
*
* \tparam BitVector A bit vector type suitable for the global indices of the basis
* \tparam Vector A coefficient vector type suitable for the global indices of the basis
* \tparam Basis Type of the global basis
* \tparam Function Type of the function to interpolate
*
* \param basis A global basis
* \param f The Dirichlet values given as a function
* \param boundaryPatch The BoundaryPath where the boundary condition should be applied
*
* This will interpolate the function `f` on the `boundaryPatch`
* return an `AffineConstraints` object constraining all corresponding DOFs.
*
* If the template default value `std::monostate` is used for the template
* parameters `BitVector` or `Vector` then the actual container types are deduced
* using `Dune::Functions::makeContainer(...)` and `Dune::Functions::makeISTLVector(...)`
* from the container descriptor provided by `basis` (only with dune-function>=2.11).
*
* \warning This function is still experimental and the interface may change.
*/
template
<
class
BitVector
=
std
::
monostate
,
class
Vector
=
std
::
monostate
,
class
Basis
,
class
Function
>
auto
makeBoundaryConstraints
(
const
Basis
&
basis
,
Function
&&
f
,
const
BoundaryPatch
<
typename
Basis
::
GridView
>&
boundaryPatch
)
{
auto
constraints
=
makeAffineConstraints
<
BitVector
,
Vector
>
(
basis
);
computeBoundaryConstraints
(
constraints
,
basis
,
f
,
boundaryPatch
);
return
constraints
;
}
}
// namespace Dune::Fufem
#endif // DUNE_FUFEM_CONSTRAINTS_BOUNDARYCONSTRAINTS_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