~njansson/dolfin/hpc

« back to all changes in this revision

Viewing changes to site-packages/dolfin/assemble.py

  • Committer: Kent-Andre Mardal
  • Date: 2008-06-07 09:10:29 UTC
  • mfrom: (2668.14.7 trunk)
  • mto: (2668.1.50 trunk)
  • mto: This revision was merged to the branch mainline in revision 2670.
  • Revision ID: kent-and@simula.no-20080607091029-ulyrrbk8kqd9xjw0
merge

Show diffs side-by-side

added added

removed removed

Lines of Context:
20
20
           "AvgMeshSize",
21
21
           "LinearPDE",
22
22
           "DirichletBC",
23
 
           "PeriodicBC"]
 
23
           "PeriodicBC",
 
24
           "compile_functions"]
24
25
 
25
26
from ffc import *
26
27
from dolfin import *
 
28
from compile_functions import compile_functions
27
29
 
28
30
# Cache for dof maps
29
31
_dof_map_cache = {}
38
40
    "Assemble form over mesh and return tensor"
39
41
 
40
42
    # Create empty list of coefficients, filled below
41
 
    coefficients_ = ArrayFunctionPtr()
 
43
    _coefficients = ArrayFunctionPtr()
42
44
 
43
45
    # Check if we need to compile the form (JIT)
44
46
    if hasattr(form, "create_cell_integral"):
47
49
        
48
50
        # Extract coefficients
49
51
        if not coefficients is None: 
50
 
            for c in coefficients:
51
 
                # Note: We could generalize this to support more objects like sympy expressions, tuples for constant vectors, etc...
 
52
            # Compile all strings as dolfin::Function
 
53
            string_expressions = []
 
54
            for c in coefficients:
 
55
                # Note: To allow tuples of floats or ints below, this logic becomes more involved...
 
56
                if isinstance(c, (tuple, str)):
 
57
                    string_expressions.append(c)
 
58
            if string_expressions:
 
59
                compiled_functions = compile_functions(string_expressions, mesh)
 
60
                compiled_functions.reverse()
 
61
            
 
62
            # Build list of coefficients
 
63
            for c in coefficients:
 
64
                # Note: We could generalize this to support more objects 
 
65
                # like sympy expressions, tuples for constant vectors, etc...
52
66
                if isinstance(c, (float, int)):
53
67
                    c = cpp_Function(mesh, float(c))
54
 
                coefficients_.push_back(c)
 
68
                elif isinstance(c, (tuple, str)):
 
69
                    c = compiled_functions.pop()
 
70
                _coefficients.push_back(c)
55
71
    else:
56
72
        # FFC form, call JIT compile
57
73
        optimize = dolfin_get("optimize form") or dolfin_get("optimize")
59
75
        
60
76
        # Extract coefficients
61
77
        for c in form_data.coefficients:
62
 
            coefficients_.push_back(c.f)
 
78
            _coefficients.push_back(c.f)
63
79
 
64
80
    # FIXME: do we need these lines now? None works fine with Assembler
65
81
    # Create dummy arguments for domains (not yet supported in Python)
82
98
        reset_tensor = True
83
99
 
84
100
    # Assemble tensor from compiled form
85
 
    cpp_assemble(tensor, compiled_form, mesh, coefficients_, dof_maps,
 
101
    cpp_assemble(tensor, compiled_form, mesh, _coefficients, dof_maps,
86
102
                 cell_domains, exterior_facet_domains, interior_facet_domains, reset_tensor)
87
103
    
88
104
    # Convert to float for scalars