This file is indexed.

/usr/share/pyshared/dolfin/fem/solving.py is in python-dolfin 1.0.0-7.

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
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
"""This module provides a small Python layer on top of the C++
VariationalProblem/Solver classes as well as the solve function."""

# Copyright (C) 2011 Anders Logg
#
# This file is part of DOLFIN.
#
# DOLFIN 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.
#
# DOLFIN 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 DOLFIN. If not, see <http://www.gnu.org/licenses/>.
#
# Modified by Marie E. Rognes, 2011.
# Modified by Johan Hake, 2011.
#
# First added:  2011-06-22
# Last changed: 2011-11-10

__all__ = ["LinearVariationalProblem",
           "LinearVariationalSolver",
           "NonlinearVariationalProblem",
           "NonlinearVariationalSolver",
           "solve"]

# Import SWIG-generated extension module (DOLFIN C++)
import dolfin.cpp as cpp

# Import UFL
import ufl

# Local imports
from dolfin.fem.form import Form
from dolfin.fem.formmanipulations import derivative

# Problem classes need special handling since they involve JIT compilation

class LinearVariationalProblem(cpp.LinearVariationalProblem):

    # Reuse C++ doc-string
    __doc__ = cpp.LinearVariationalProblem.__doc__

    def __init__(self, a, L, u, bcs=None,
                 form_compiler_parameters=None):
        """
        Create linear variational problem a(u, v) = L(v).

        An optional argument bcs may be passed to specify boundary
        conditions.

        Another optional argument form_compiler_parameters may be
        specified to pass parameters to the form compiler.
        """

        # Extract and check arguments
        u = _extract_u(u)
        bcs = _extract_bcs(bcs)

        # Store input UFL forms and solution Function
        self.a_ufl = a
        self.L_ufl = L
        self.u_ufl = u

        # Store form compiler parameters
        form_compiler_parameters = form_compiler_parameters or {}
        self.form_compiler_parameters = form_compiler_parameters

        # Wrap forms (and check if linear form L is empty)
        if L.integrals() == ():
            L = cpp.Form(1, 0)
        else:
            L = Form(L, form_compiler_parameters=form_compiler_parameters)
        a = Form(a, form_compiler_parameters=form_compiler_parameters)

        # Initialize C++ base class
        cpp.LinearVariationalProblem.__init__(self, a, L, u, bcs)

class NonlinearVariationalProblem(cpp.NonlinearVariationalProblem):

    # Reuse C++ doc-string
    __doc__ = cpp.NonlinearVariationalProblem.__doc__

    def __init__(self, F, u, bcs=None, J=None,
                 form_compiler_parameters=None):
        """
        Create nonlinear variational problem F(u; v) = 0.

        Optional arguments bcs and J may be passed to specify boundary
        conditions and the Jacobian J = dF/du.

        Another optional argument form_compiler_parameters may be
        specified to pass parameters to the form compiler.
        """

        # Extract and check arguments
        u = _extract_u(u)
        bcs = _extract_bcs(bcs)

        # Store input UFL forms and solution Function
        self.F_ufl = F
        self.J_ufl = J
        self.u_ufl = u

        # Store form compiler parameters
        form_compiler_parameters = form_compiler_parameters or {}
        self.form_compiler_parameters = form_compiler_parameters

        # Wrap forms
        F = Form(F, form_compiler_parameters=form_compiler_parameters)
        if J is not None:
            J = Form(J, form_compiler_parameters=form_compiler_parameters)

        # Initialize C++ base class
        if J is None:
            cpp.NonlinearVariationalProblem.__init__(self, F, u, bcs)
        else:
            cpp.NonlinearVariationalProblem.__init__(self, F, u, bcs, J)

# Solver classes are imported directly
from dolfin.cpp import LinearVariationalSolver, NonlinearVariationalSolver
from dolfin.fem.adaptivesolving import AdaptiveLinearVariationalSolver
from dolfin.fem.adaptivesolving import AdaptiveNonlinearVariationalSolver

# Solve function handles both linear systems and variational problems

def solve(*args, **kwargs):
    """Solve linear system Ax = b or variational problem a == L or F == 0.

    The DOLFIN solve() function can be used to solve either linear
    systems or variational problems. The following list explains the
    various ways in which the solve() function can be used.

    *1. Solving linear systems*

    A linear system Ax = b may be solved by calling solve(A, x, b),
    where A is a matrix and x and b are vectors. Optional arguments
    may be passed to specify the solver method and preconditioner.
    Some examples are given below:

    .. code-block:: python

        solve(A, x, b)
        solve(A, x, b, "lu")
        solve(A, x, b, "gmres", "ilu")
        solve(A, x, b, "cg", "hypre_amg")

    Possible values for the solver method and preconditioner depend
    on which linear algebra backend is used and how that has been
    configured.

    To list all available LU methods, run the following command:

    .. code-block:: python

        list_lu_methods()

    To list all available Krylov methods, run the following command:

    .. code-block:: python

        list_krylov_methods()

    To list all available preconditioners, run the following command:

    .. code-block:: python

        list_preconditioners()

    To list all available solver methods, including LU methods, Krylov
    methods and, possibly, other methods, run the following command:

    .. code-block:: python

        list_solver_methods()

    *2. Solving linear variational problems*

    A linear variational problem a(u, v) = L(v) for all v may be
    solved by calling solve(a == L, u, ...), where a is a bilinear
    form, L is a linear form, u is a Function (the solution). Optional
    arguments may be supplied to specify boundary conditions or solver
    parameters. Some examples are given below:

    .. code-block:: python

        solve(a == L, u)
        solve(a == L, u, bcs=bc)
        solve(a == L, u, bcs=[bc1, bc2])

        solve(a == L, u, bcs=bcs,
              solver_parameters={"linear_solver": "lu"},
              form_compiler_parameters={"optimize": True})

    For available choices for the 'solver_parameters' kwarg, look at:

        .. code-block:: python

        info(LinearVariationalSolver.default_parameters(), 1)

    *3. Solving nonlinear variational problems*

    A nonlinear variational problem F(u; v) = 0 for all v may be
    solved by calling solve(F == 0, u, ...), where the residual F is a
    linear form (linear in the test function v but possibly nonlinear
    in the unknown u) and u is a Function (the solution). Optional
    arguments may be supplied to specify boundary conditions, the
    Jacobian form or solver parameters. If the Jacobian is not
    supplied, it will be computed by automatic differentiation of the
    residual form. Some examples are given below:

    .. code-block:: python

        solve(F == 0, u)
        solve(F == 0, u, bcs=bc)
        solve(F == 0, u, bcs=[bc1, bc2])

        solve(F == 0, u, bcs, J=J,
              solver_parameters={"linear_solver": "lu"},
              form_compiler_parameters={"optimize": True})


    For available choices for the 'solver_parameters' kwarg, look at:

        .. code-block:: python

        info(NonlinearVariationalSolver.default_parameters(), 1)

    *4. Solving linear/nonlinear variational problems adaptively

    FIXME: Missing documentation for adaptive solve, add here.


    """

    assert(len(args) > 0)

    # Call adaptive solve if we get a tolerance
    if "tol" in kwargs:
        _solve_varproblem_adaptive(*args, **kwargs)

    # Call variational problem solver if we get an equation (but not a tolerance)
    elif isinstance(args[0], ufl.classes.Equation):
        _solve_varproblem(*args, **kwargs)

    # Default case, just call the wrapped C++ solve function
    else:
        return cpp.solve(*args)

def _solve_varproblem(*args, **kwargs):
    "Solve variational problem a == L or F == 0"

    # Extract arguments
    eq, u, bcs, J, tol, M, form_compiler_parameters, solver_parameters \
        = _extract_args(*args, **kwargs)

    # Solve linear variational problem
    if isinstance(eq.lhs, ufl.Form) and isinstance(eq.rhs, ufl.Form):

        # Create problem
        problem = LinearVariationalProblem(eq.lhs, eq.rhs, u, bcs,
                    form_compiler_parameters=form_compiler_parameters)

        # Create solver and call solve
        solver = LinearVariationalSolver(problem)
        solver.parameters.update(solver_parameters)
        solver.solve()

    # Solve nonlinear variational problem
    else:

        # Create Jacobian if missing
        if J is None:
            cpp.info("No Jacobian form specified for nonlinear variational problem.")
            cpp.info("Differentiating residual form F to obtain Jacobian J = F'.")
            F = eq.lhs
            J = derivative(F, u)

        # Create problem
        problem = NonlinearVariationalProblem(eq.lhs, u, bcs, J,
                    form_compiler_parameters=form_compiler_parameters)

        # Create solver and call solve
        solver = NonlinearVariationalSolver(problem)
        solver.parameters.update(solver_parameters)
        solver.solve()

def _solve_varproblem_adaptive(*args, **kwargs):
    "Solve variational problem a == L or F == 0 adaptively"

    # Extract arguments
    eq, u, bcs, J, tol, M, form_compiler_parameters, solver_parameters \
        = _extract_args(*args, **kwargs)

    # Check that we received the goal functional
    if M is None:
        cpp.dolfin_error("solving.py",
                         "solve variational problem adaptively",
                         "Missing goal functional")

    # Solve linear variational problem
    if isinstance(eq.lhs, ufl.Form) and isinstance(eq.rhs, ufl.Form):

        # Create problem
        problem = LinearVariationalProblem(eq.lhs, eq.rhs, u, bcs,
                    form_compiler_parameters=form_compiler_parameters)

        # Create solver and call solve
        solver = AdaptiveLinearVariationalSolver(problem)
        solver.parameters.update(solver_parameters)
        solver.solve(tol, M)

    # Solve nonlinear variational problem
    else:

        # Create Jacobian if missing
        if J is None:
            cpp.info("No Jacobian form specified for nonlinear variational problem.")
            cpp.info("Differentiating residual form F to obtain Jacobian J = F'.")
            F = eq.lhs
            J = derivative(F, u)

        # Create problem
        problem = NonlinearVariationalProblem(eq.lhs, u, bcs, J,
                    form_compiler_parameters=form_compiler_parameters)

        # Create solver and call solve
        solver = AdaptiveNonlinearVariationalSolver(problem)
        solver.parameters.update(solver_parameters)
        solver.solve(tol, M)

def _extract_args(*args, **kwargs):
    "Common extraction of arguments for _solve_varproblem[_adaptive]"

    # Check for use of valid kwargs
    valid_kwargs = ["bcs", "J", "tol", "M",
                    "form_compiler_parameters", "solver_parameters"]
    for kwarg in kwargs.iterkeys():
        if not kwarg in valid_kwargs:
            cpp.dolfin_error("solving.py",
                             "solve variational problem",
                             "Illegal keyword argument \"%s\"; valid keywords are %s" % \
                                 (kwarg,
                                  ", ".join("\"%s\"" % kwarg for kwarg in valid_kwargs)))

    # Extract equation
    if not len(args) >= 2:
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Missing arguments, expecting solve(lhs == rhs, "\
                         "u, bcs=bcs), where bcs is optional")
    if len(args) > 3:
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Too many arguments, expecting solve(lhs == rhs, "\
                         "u, bcs=bcs), where bcs is optional")

    # Extract equation
    eq = _extract_eq(args[0])

    # Extract solution function
    u = _extract_u(args[1])

    # Extract boundary conditions
    if len(args) > 2:
        bcs = _extract_bcs(args[2])
    elif "bcs" in kwargs:
        bcs = _extract_bcs(kwargs["bcs"])
    else:
        bcs = []

    # Extract Jacobian
    J = kwargs.get("J", None)
    if J is not None and not isinstance(J, ufl.Form):
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Expecting Jacobian J to be a UFL Form")

    # Extract tolerance
    tol = kwargs.get("tol", None)
    if tol is not None and not (isinstance(tol, (float, int)) and tol >= 0.0):
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Expecting tolerance tol to be a non-negative number")

    # Extract functional
    M = kwargs.get("M", None)
    if M is not None and not isinstance(M, ufl.Form):
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Expecting goal functional M to be a UFL Form")

    # Extract parameters
    form_compiler_parameters = kwargs.get("form_compiler_parameters", {})
    solver_parameters = kwargs.get("solver_parameters", {})

    return eq, u, bcs, J, tol, M, form_compiler_parameters, solver_parameters

def _extract_eq(eq):
    "Extract and check argument eq"
    if not isinstance(eq, ufl.classes.Equation):
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Expecting first argument to be an Equation")
    return eq

def _extract_u(u):
    "Extract and check argument u"
    if not isinstance(u, cpp.Function):
        cpp.dolfin_error("solving.py",
                         "solve variational problem",
                         "Expecting second argument to be a Function")
    return u

def _extract_bcs(bcs):
    "Extract and check argument bcs"
    if bcs is None:
        bcs = []
    elif not isinstance(bcs, (list, tuple)):
        bcs = [bcs]
    for bc in bcs:
        if not isinstance(bc, cpp.BoundaryCondition):
            cpp.dolfin_error("solving.py",
                             "solve variational problem",
                             "Unable to extract boundary condition arguments")
    return bcs