/usr/lib/python2.7/dist-packages/ffc/errorcontrol/errorcontrolgenerators.py is in python-ffc 1.3.0-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 | """
This module provides an abstract ErrorControlGenerator class for
generating forms required for goal-oriented error control and a
realization of this: UFLErrorControlGenerator for handling pure UFL
forms.
"""
# Copyright (C) 2010 Marie E. Rognes
#
# This file is part of FFC.
#
# FFC is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# FFC is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with FFC. If not, see <http://www.gnu.org/licenses/>.
#
# Last changed: 2011-06-29
from ufl import inner, dx, ds, dS, avg, replace, action
from ufl.algorithms.analysis import extract_arguments
__all__ = ["ErrorControlGenerator", "UFLErrorControlGenerator"]
class ErrorControlGenerator:
def __init__(self, module, F, M, u):
"""
*Arguments*
module (Python module)
The module to use for specific form manipulations
(typically ufl or dolfin)
F (tuple or Form)
tuple of (bilinear, linear) forms or linear form
M (Form)
functional or linear form
u (Coefficient)
The coefficient considered as the unknown.
"""
# Store module
self.module = module
# Store solution Coefficient/Function
self.u = u
# Extract the lhs (bilinear form), rhs (linear form), goal
# (functional), weak residual (linear form)
linear_case = (isinstance(F, (tuple, list)) and len(F) == 2)
if linear_case:
self.lhs, self.rhs = F
try:
self.goal = action(M, u)
except:
self.goal = M # Allow functionals as input as well
self.weak_residual = self.rhs - action(self.lhs, u)
else:
self.lhs = self.module.derivative(F, u)
self.rhs = F
self.goal = M
self.weak_residual = - F
# At least check that final forms have correct rank
assert(len(extract_arguments(self.lhs)) == 2)
assert(len(extract_arguments(self.rhs)) == 1)
assert(len(extract_arguments(self.goal)) == 0)
assert(len(extract_arguments(self.weak_residual)) == 1)
# Store map from identifiers to names for forms and generated
# coefficients
self.ec_names = {}
# Use predefined names for the forms in the primal problem
self.ec_names[id(self.lhs)] = "lhs"
self.ec_names[id(self.rhs)] = "rhs"
self.ec_names[id(self.goal)] = "goal"
# Initialize other required data
self.initialize_data()
def initialize_data(self):
"""
Initialize specific data
"""
msg = """ErrorControlGenerator acts as an abstract
class. Subclasses must overload the initialize_data() method
and provide a certain set of variables. See
UFLErrorControlGenerator for an example."""
raise NotImplementedError, msg
def generate_all_error_control_forms(self):
"""
Utility function for generating all (8) forms required for
error control in addition to the primal forms
"""
# Generate dual forms
(a_star, L_star) = self.dual_forms()
# Generate forms for computing strong cell residual
(a_R_T, L_R_T) = self.cell_residual()
# Generate forms for computing strong facet residuals
(a_R_dT, L_R_dT) = self.facet_residual()
# Generate form for computing error estimate
eta_h = self.error_estimate()
# Generate form for computing error indicators
eta_T = self.error_indicators()
# Return all generated forms in CERTAIN order matching
# constructor of dolfin/adaptivity/ErrorControl.h
return (a_star, L_star, eta_h, a_R_T, L_R_T, a_R_dT, L_R_dT, eta_T)
def primal_forms(self):
"""
Return primal forms in order (bilinear, linear, functional)
"""
return self.lhs, self.rhs, self.goal
def dual_forms(self):
"""
Generate and return (bilinear, linear) forms defining linear
dual variational problem
"""
a_star = self.module.adjoint(self.lhs)
L_star = self.module.derivative(self.goal, self.u)
return (a_star, L_star)
def cell_residual(self):
"""
Generate and return (bilinear, linear) forms defining linear
variational problem for the strong cell residual
"""
# Define trial and test functions for the cell residuals on
# discontinuous version of primal trial space
R_T = self.module.TrialFunction(self._dV)
v = self.module.TestFunction(self._dV)
# Extract original test function in the weak residual
v_h = extract_arguments(self.weak_residual)[0]
# Define forms defining linear variational problem for cell
# residual
v_T = self._b_T*v
a_R_T = inner(v_T, R_T)*dx()
L_R_T = replace(self.weak_residual, {v_h: v_T})
return (a_R_T, L_R_T)
def facet_residual(self):
"""
Generate and return (bilinear, linear) forms defining linear
variational problem for the strong facet residual(s)
"""
# Define trial and test functions for the facet residuals on
# discontinuous version of primal trial space
R_e = self.module.TrialFunction(self._dV)
v = self.module.TestFunction(self._dV)
# Extract original test function in the weak residual
v_h = extract_arguments(self.weak_residual)[0]
# Define forms defining linear variational problem for facet
# residual
v_e = self._b_e*v
a_R_dT = ((inner(v_e('+'), R_e('+')) + inner(v_e('-'), R_e('-')))*dS()
+ inner(v_e, R_e)*ds())
L_R_dT = (replace(self.weak_residual, {v_h: v_e})
- inner(v_e, self._R_T)*dx())
return (a_R_dT, L_R_dT)
def error_estimate(self):
"""
Generate and return functional defining error estimate
"""
# Error estimate defined as r(Ez_h):
eta_h = action(self.weak_residual, self._Ez_h)
return eta_h
def error_indicators(self):
"""
Generate and return linear form defining error indicators
"""
# Extract these to increase readability
R_T = self._R_T
R_dT = self._R_dT
z = self._Ez_h
z_h = self._z_h
# Define linear form for computing error indicators
v = self.module.TestFunction(self._DG0)
eta_T = (v*inner(R_T, z - z_h)*dx()
+ avg(v)*(inner(R_dT('+'), (z - z_h)('+'))
+ inner(R_dT('-'), (z - z_h)('-')))*dS()
+ v*inner(R_dT, z - z_h)*ds())
return eta_T
class UFLErrorControlGenerator(ErrorControlGenerator):
"""
This class provides a realization of ErrorControlGenerator for use
with pure UFL forms
"""
def __init__(self, F, M, u):
"""
*Arguments*
F (tuple or Form)
tuple of (bilinear, linear) forms or linear form
M (Form)
functional or linear form
u (Coefficient)
The coefficient considered as the unknown.
"""
ErrorControlGenerator.__init__(self, __import__("ufl"), F, M, u)
def initialize_data(self):
"""
Extract required objects for defining error control
forms. This will be stored, reused and in particular named.
"""
# Developer's note: The UFL-FFC-DOLFIN--PyDOLFIN toolchain for
# error control is quite fine-tuned. In particular, the order
# of coefficients in forms is (and almost must be) used for
# their assignment. This means that the order in which these
# coefficients are defined matters and should be considered
# fixed.
from ufl import FiniteElement, Coefficient
from ufl.algorithms.elementtransformations import tear, increase_order
# Primal trial element space
self._V = self.u.element()
# Primal test space == Dual trial space
Vhat = extract_arguments(self.weak_residual)[0].element()
# Discontinuous version of primal trial element space
self._dV = tear(self._V)
# Extract cell and geometric dimension
cell = self._V.cell()
gdim = cell.geometric_dimension()
# Coefficient representing improved dual
E = increase_order(Vhat)
self._Ez_h = Coefficient(E)
self.ec_names[id(self._Ez_h)] = "__improved_dual"
# Coefficient representing cell bubble function
B = FiniteElement("B", cell, gdim + 1)
self._b_T = Coefficient(B)
self.ec_names[id(self._b_T)] = "__cell_bubble"
# Coefficient representing strong cell residual
self._R_T = Coefficient(self._dV)
self.ec_names[id(self._R_T)] = "__cell_residual"
# Coefficient representing cell cone function
C = FiniteElement("DG", cell, gdim)
self._b_e = Coefficient(C)
self.ec_names[id(self._b_e)] = "__cell_cone"
# Coefficient representing strong facet residual
self._R_dT = Coefficient(self._dV)
self.ec_names[id(self._R_dT)] = "__facet_residual"
# Define discrete dual on primal test space
self._z_h = Coefficient(Vhat)
self.ec_names[id(self._z_h)] = "__discrete_dual_solution"
# Piecewise constants for assembling indicators
self._DG0 = FiniteElement("DG", cell, 0)
|