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
Iterations
Wiki
Requirements
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Locked files
Build
Pipelines
Jobs
Pipeline schedules
Test cases
Artifacts
Deploy
Releases
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
This is an archived project. Repository and other project resources are read-only.
Show more breadcrumbs
agnumpde
dune-fufem
Commits
6e9f47b6
Commit
6e9f47b6
authored
10 years ago
by
Max Kahnt
Browse files
Options
Downloads
Patches
Plain Diff
Remove necessity to specify fallback QuadRuleKey in SubgridL2FunctionalAssembler.
Provide two separate constructors instead.
parent
6c5c5563
Branches
Branches containing commit
No related tags found
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dune/fufem/assemblers/localassemblers/subgridl2functionalassembler.hh
+32
-13
32 additions, 13 deletions
...ssemblers/localassemblers/subgridl2functionalassembler.hh
with
32 additions
and
13 deletions
dune/fufem/assemblers/localassemblers/subgridl2functionalassembler.hh
+
32
−
13
View file @
6e9f47b6
...
@@ -50,22 +50,41 @@ class SubgridL2FunctionalAssembler :
...
@@ -50,22 +50,41 @@ class SubgridL2FunctionalAssembler :
* It can assemble functionals on the subgrid given by grid
* It can assemble functionals on the subgrid given by grid
* functions on the underlying hostgrid exactly.
* functions on the underlying hostgrid exactly.
*
*
* The QuadratureRuleKey given here does only specify what is
* \param f the basisgrid function representing the functional
* needed to integrate f if f is not a BasisGridFunction on
* the hostgrid. Otherwise a quadrature rule is determined
* automatically.
*
* \param f the (hostgrid) function representing the functional
* \param grid the subgrid (!)
* \param grid the subgrid (!)
* \param fQuadKey A fallback QuadratureRuleKey that specifies how to integrate f if its not a BasisGridFunction on the hostgrid
*/
*/
SubgridL2FunctionalAssembler
(
const
HostGridFunction
&
f
,
const
GridType
&
grid
,
const
QuadratureRuleKey
&
fQuadKey
)
:
template
<
typename
B
,
typename
C
,
std
::
enable_if
<
std
::
is_base_of
<
HostGridFunction
,
BasisGridFunction
<
B
,
C
>
>::
value_type
>>
functionQuadKey_
(
fQuadKey
),
SubgridL2FunctionalAssembler
(
const
BasisGridFunction
<
B
,
C
>&
f
,
const
GridType
&
grid
)
:
functionQuadKey_
(
0
),
grid_
(
grid
),
grid_
(
grid
),
f_
(
f
),
f_
(
f
),
hostBasisGridFunctionInfo_
(
dynamic_cast
<
const
HostBasisGridFunctionInfo
*>
(
&
f
)
)
hostBasisGridFunctionInfo_
(
&
f
)
{}
{}
/** \brief constructor
*
* Creates a local functional assembler for an L2-functional.
* It can assemble functionals on the subgrid given by grid
* functions on the underlying hostgrid exactly.
*
* The QuadratureRuleKey given here does specify what is
* needed to integrate f. It can be used to override the
* automatically determined QuadratureRule if f is a BasisGridFunction.
*
* \param f the hostgrid function representing the functional
* \param grid the subgrid (!)
* \param fQuadKey A QuadratureRuleKey that specifies how to integrate f
*/
SubgridL2FunctionalAssembler
(
const
HostGridFunction
&
f
,
const
GridType
&
grid
,
const
QuadratureRuleKey
&
fQuadKey
)
:
functionQuadKey_
(
&
fQuadKey
),
grid_
(
grid
),
f_
(
f
),
hostBasisGridFunctionInfo_
(
0
)
{
if
(
dynamic_cast
<
const
HostBasisGridFunctionInfo
*>
(
&
f
))
Dune
::
dwarn
<<
"Overriding HostBasis Quadrature rule."
<<
std
::
endl
;
}
/** \copydoc L2FunctionalAssembler::assemble
/** \copydoc L2FunctionalAssembler::assemble
*/
*/
void
assemble
(
const
Element
&
element
,
LocalVector
&
localVector
,
const
TrialLocalFE
&
tFE
)
const
void
assemble
(
const
Element
&
element
,
LocalVector
&
localVector
,
const
TrialLocalFE
&
tFE
)
const
...
@@ -90,7 +109,7 @@ class SubgridL2FunctionalAssembler :
...
@@ -90,7 +109,7 @@ class SubgridL2FunctionalAssembler :
if
(
hostBasisGridFunctionInfo_
)
if
(
hostBasisGridFunctionInfo_
)
quadKey
=
quadKey
.
product
(
hostBasisGridFunctionInfo_
->
quadratureRuleKey
(
*
baseHostElement
));
quadKey
=
quadKey
.
product
(
hostBasisGridFunctionInfo_
->
quadratureRuleKey
(
*
baseHostElement
));
else
else
quadKey
=
quadKey
.
product
(
functionQuadKey_
);
quadKey
=
quadKey
.
product
(
*
functionQuadKey_
);
const
Dune
::
template
QuadratureRule
<
double
,
dim
>
&
quad
=
QuadratureRuleCache
<
double
,
dim
>::
rule
(
quadKey
);
const
Dune
::
template
QuadratureRule
<
double
,
dim
>
&
quad
=
QuadratureRuleCache
<
double
,
dim
>::
rule
(
quadKey
);
// loop over quadrature points
// loop over quadrature points
...
@@ -136,7 +155,7 @@ class SubgridL2FunctionalAssembler :
...
@@ -136,7 +155,7 @@ class SubgridL2FunctionalAssembler :
if
(
hostBasisGridFunctionInfo_
)
if
(
hostBasisGridFunctionInfo_
)
quadKey
=
quadKey
.
product
(
hostBasisGridFunctionInfo_
->
quadratureRuleKey
(
hostElement
));
quadKey
=
quadKey
.
product
(
hostBasisGridFunctionInfo_
->
quadratureRuleKey
(
hostElement
));
else
else
quadKey
=
quadKey
.
product
(
functionQuadKey_
);
quadKey
=
quadKey
.
product
(
*
functionQuadKey_
);
const
Dune
::
template
QuadratureRule
<
double
,
dim
>
&
quad
=
QuadratureRuleCache
<
double
,
dim
>::
rule
(
quadKey
);
const
Dune
::
template
QuadratureRule
<
double
,
dim
>
&
quad
=
QuadratureRuleCache
<
double
,
dim
>::
rule
(
quadKey
);
...
@@ -171,7 +190,7 @@ class SubgridL2FunctionalAssembler :
...
@@ -171,7 +190,7 @@ class SubgridL2FunctionalAssembler :
}
}
private
:
private
:
const
QuadratureRuleKey
functionQuadKey_
;
const
QuadratureRuleKey
*
functionQuadKey_
;
const
GridType
&
grid_
;
const
GridType
&
grid_
;
const
HostGridFunction
&
f_
;
const
HostGridFunction
&
f_
;
const
HostBasisGridFunctionInfo
*
hostBasisGridFunctionInfo_
;
const
HostBasisGridFunctionInfo
*
hostBasisGridFunctionInfo_
;
...
...
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