Coverage for python/src/dolfinx_mpc/assemble_vector.py: 100%

58 statements  

« prev     ^ index     » next       coverage.py v7.10.7, created at 2025-09-29 19:56 +0000

1# Copyright (C) 2020 Jørgen S. Dokken 

2# 

3# This file is part of DOLFINX_MPC 

4# 

5# SPDX-License-Identifier: MIT 

6from __future__ import annotations 

7 

8import contextlib 

9from typing import Iterable, Optional, Sequence, Union 

10 

11from petsc4py import PETSc as _PETSc 

12 

13import dolfinx.cpp as _cpp 

14import dolfinx.fem as _fem 

15import numpy 

16import ufl 

17from dolfinx import default_scalar_type 

18from dolfinx.common import Timer 

19from dolfinx.la.petsc import create_vector 

20 

21import dolfinx_mpc.cpp 

22 

23from .multipointconstraint import MultiPointConstraint, _float_classes 

24 

25 

26def apply_lifting( 

27 b: _PETSc.Vec, # type: ignore 

28 form: Union[Iterable[Sequence[_fem.Form]], Iterable[_fem.Form]], # type: ignore 

29 bcs: Union[Sequence[_fem.DirichletBC], Sequence[Sequence[_fem.DirichletBC]]], 

30 constraint: Union[MultiPointConstraint, Sequence[MultiPointConstraint]], 

31 x0: Optional[Sequence[_PETSc.Vec]] = None, # type: ignore 

32 scale: _float_classes = default_scalar_type(1.0), 

33): # type: ignore 

34 """ 

35 Apply lifting to vector b, i.e. 

36 :math:`b = b - scale \\cdot K^T (A_j (g_j - x0_j))` 

37 

38 Args: 

39 b: PETSc vector to assemble into 

40 form: The linear form 

41 bcs: List of Dirichlet boundary conditions 

42 constraint: The multi point constraint 

43 x0: List of vectors 

44 scale: Scaling for lifting 

45 """ 

46 t = Timer("~MPC: Apply lifting (C++)") 

47 if isinstance(scale, numpy.generic): # nanobind conversion of numpy dtypes to general Python types 

48 scale = scale.item() # type: ignore 

49 

50 if b.getType() == "nest": 

51 try: 

52 bcs = _fem.bcs_by_block(_fem.extract_function_spaces(form, 1), bcs) # type: ignore 

53 except AttributeError: 

54 pass 

55 x0 = [] if x0 is None else x0.getNestSubVecs() # type: ignore 

56 assert isinstance(form, Sequence) and isinstance(constraint, Sequence) 

57 for b_sub, a_sub, mpc_i in zip(b.getNestSubVecs(), form, constraint): 

58 _a = [None if form is None else form._cpp_object for form in a_sub] # type:ignore 

59 _bcs = [[bc._cpp_object for bc in bcs0] for bcs0 in bcs] # type: ignore 

60 dolfinx_mpc.cpp.mpc.apply_lifting(b_sub.array_w, _a, _bcs, x0, scale, mpc_i._cpp_object) 

61 else: 

62 with contextlib.ExitStack() as stack: 

63 if x0 is None: 

64 x0 = [] 

65 else: 

66 x0 = [stack.enter_context(x.localForm()) for x in x0] 

67 x0_r = [x.array_r for x in x0] 

68 b_local = stack.enter_context(b.localForm()) 

69 _forms = [f._cpp_object for f in form] # type: ignore 

70 _bcs = [[bc._cpp_object for bc in bcs0] for bcs0 in bcs] # type: ignore 

71 assert isinstance(constraint, MultiPointConstraint) 

72 dolfinx_mpc.cpp.mpc.apply_lifting(b_local.array_w, _forms, _bcs, x0_r, scale, constraint._cpp_object) 

73 t.stop() 

74 

75 

76def assemble_vector( 

77 form: ufl.form.Form, 

78 constraint: MultiPointConstraint, 

79 b: Optional[_PETSc.Vec] = None, # type: ignore 

80) -> _PETSc.Vec: # type: ignore 

81 """ 

82 Assemble a linear form into vector `b` with corresponding multi point constraint 

83 

84 Args: 

85 form: The linear form 

86 constraint: The multi point constraint 

87 b: PETSc vector to assemble 

88 

89 Returns: 

90 The vector with the assembled linear form (`b` if supplied) 

91 """ 

92 

93 if b is None: 

94 b = create_vector([(constraint.function_space.dofmap.index_map, constraint.function_space.dofmap.index_map_bs)]) 

95 t = Timer("~MPC: Assemble vector (C++)") 

96 with b.localForm() as b_local: 

97 b_local.set(0.0) 

98 dolfinx_mpc.cpp.mpc.assemble_vector(b_local.array_w, form._cpp_object, constraint._cpp_object) 

99 t.stop() 

100 return b 

101 

102 

103def create_vector_nest(L: Sequence[_fem.Form], constraints: Sequence[MultiPointConstraint]) -> _PETSc.Vec: # type: ignore 

104 """ 

105 Create a PETSc vector of type "nest" appropriate for the provided multi 

106 point constraints 

107 

108 Args: 

109 L: A sequence of linear forms 

110 constraints: An ordered list of multi point constraints 

111 

112 Returns: 

113 PETSc.Vec: A PETSc vector of type "nest" #type: ignore 

114 """ 

115 assert len(constraints) == len(L) 

116 

117 maps = [ 

118 (constraint.function_space.dofmap.index_map, constraint.function_space.dofmap.index_map_bs) 

119 for constraint in constraints 

120 ] 

121 return _cpp.fem.petsc.create_vector_nest(maps) 

122 

123 

124def assemble_vector_nest( 

125 b: _PETSc.Vec, # type: ignore 

126 L: Sequence[_fem.Form], 

127 constraints: Sequence[MultiPointConstraint], 

128): 

129 """ 

130 Assemble a linear form into a PETSc vector of type "nest" 

131 

132 Args: 

133 b: A PETSc vector of type "nest" 

134 L: A sequence of linear forms 

135 constraints: An ordered list of multi point constraints 

136 """ 

137 assert len(constraints) == len(L) 

138 assert b.getType() == "nest" 

139 

140 b_sub_vecs = b.getNestSubVecs() 

141 for i, L_row in enumerate(L): 

142 assemble_vector(L_row, constraints[i], b=b_sub_vecs[i])