gel.fix_dofs_overloaded
Overriden implementation of control variable restriction function for automatic differentiation in pyadjoint.
Overrides gel.fix_dofs.fix_dofs
with a Block object.
Credit https://fenicsproject.discourse.group/t/specific-dof-of-a-function-as-control-variable/3300/3
1"""Overriden implementation of control variable restriction function for 2automatic differentiation in pyadjoint. 3 4Overrides `gel.fix_dofs.fix_dofs` with a Block object. 5 6Credit https://fenicsproject.discourse.group/t/specific-dof-of-a-function-as-control-variable/3300/3""" 7from fenics import * 8from fenics_adjoint import * 9from pyadjoint import Block 10from pyadjoint.overloaded_function import overload_function 11 12from .fix_dofs import fix_dofs 13 14 15_backend_fix_dofs = fix_dofs 16 17class FixedDofsBlock(Block): 18 """Tape Block for `gel.fix_dofs.fix_dofs` 19 20 Describes how to recompute and how to compute the adjoint. 21 """ 22 23 def __init__(self, func, fixed_indexes, fixed_values, **kwargs): 24 """Take in arguments to the function and kwargs for adjoint.""" 25 super(FixedDofsBlock, self).__init__() 26 self.kwargs = kwargs 27 self.add_dependency(func) 28 29 # variables that do not appear in the dependencies but are still used 30 # for the computation 31 self.fixed_indexes = fixed_indexes 32 self.fixed_values = fixed_values 33 34 def __str__(self): 35 return 'FixedDofsBlock' 36 37 def evaluate_adj_component(self, 38 inputs, 39 adj_inputs, 40 block_variable, 41 idx, 42 prepared=None 43 ): 44 """Zero-out the components of fixed DoFs for derivative. 45 46 * `inputs` is the list of inputs to `gel.fix_dofs.fix_dofs`, ie. 47 index 0 has the function. 48 * `adj_inputs` is a list with the same number of items as 49 `inputs`, has the values of the adjoint input 50 * Other inputs unused. 51 52 Returns: the resulting product accounting for fixed DoFs 53 """ 54 adj_input = adj_inputs[0] 55 x = inputs[0].vector() 56 57 output = adj_input 58 # the derivative of the fixed dofs is null 59 output[self.fixed_indexes] = 0.0 60 return output 61 62 def recompute_component(self, inputs, block_variable, idx, prepared): 63 """Re-fixes the DoFs 64 65 * `inputs` is the list of inputs to `fix_dofs` originally 66 encountered 67 * Rest of inputs unused 68 69 Returns: the DoF-fixed function. 70 """ 71 return _backend_fix_dofs( 72 inputs[0], 73 self.fixed_indexes, 74 self.fixed_values 75 ) 76 77fix_dofs = overload_function(fix_dofs, FixedDofsBlock) 78"""Forces DoF values of a FE function to specific values. 79 80Is equipped for automatic differentiation. 81 82* `func` is a FEniCS function in a space with DoFs 83* `fixed_indexes` is an iterable of integer local DoF indices 84* `fixed_values` is an iterable of the same size as `fixed_indexes` 85with the float values to set the DoFs to. 86"""
class
FixedDofsBlock(pyadjoint.block.Block):
18class FixedDofsBlock(Block): 19 """Tape Block for `gel.fix_dofs.fix_dofs` 20 21 Describes how to recompute and how to compute the adjoint. 22 """ 23 24 def __init__(self, func, fixed_indexes, fixed_values, **kwargs): 25 """Take in arguments to the function and kwargs for adjoint.""" 26 super(FixedDofsBlock, self).__init__() 27 self.kwargs = kwargs 28 self.add_dependency(func) 29 30 # variables that do not appear in the dependencies but are still used 31 # for the computation 32 self.fixed_indexes = fixed_indexes 33 self.fixed_values = fixed_values 34 35 def __str__(self): 36 return 'FixedDofsBlock' 37 38 def evaluate_adj_component(self, 39 inputs, 40 adj_inputs, 41 block_variable, 42 idx, 43 prepared=None 44 ): 45 """Zero-out the components of fixed DoFs for derivative. 46 47 * `inputs` is the list of inputs to `gel.fix_dofs.fix_dofs`, ie. 48 index 0 has the function. 49 * `adj_inputs` is a list with the same number of items as 50 `inputs`, has the values of the adjoint input 51 * Other inputs unused. 52 53 Returns: the resulting product accounting for fixed DoFs 54 """ 55 adj_input = adj_inputs[0] 56 x = inputs[0].vector() 57 58 output = adj_input 59 # the derivative of the fixed dofs is null 60 output[self.fixed_indexes] = 0.0 61 return output 62 63 def recompute_component(self, inputs, block_variable, idx, prepared): 64 """Re-fixes the DoFs 65 66 * `inputs` is the list of inputs to `fix_dofs` originally 67 encountered 68 * Rest of inputs unused 69 70 Returns: the DoF-fixed function. 71 """ 72 return _backend_fix_dofs( 73 inputs[0], 74 self.fixed_indexes, 75 self.fixed_values 76 )
Tape Block for gel.fix_dofs.fix_dofs
Describes how to recompute and how to compute the adjoint.
FixedDofsBlock(func, fixed_indexes, fixed_values, **kwargs)
24 def __init__(self, func, fixed_indexes, fixed_values, **kwargs): 25 """Take in arguments to the function and kwargs for adjoint.""" 26 super(FixedDofsBlock, self).__init__() 27 self.kwargs = kwargs 28 self.add_dependency(func) 29 30 # variables that do not appear in the dependencies but are still used 31 # for the computation 32 self.fixed_indexes = fixed_indexes 33 self.fixed_values = fixed_values
Take in arguments to the function and kwargs for adjoint.
def
evaluate_adj_component(self, inputs, adj_inputs, block_variable, idx, prepared=None):
38 def evaluate_adj_component(self, 39 inputs, 40 adj_inputs, 41 block_variable, 42 idx, 43 prepared=None 44 ): 45 """Zero-out the components of fixed DoFs for derivative. 46 47 * `inputs` is the list of inputs to `gel.fix_dofs.fix_dofs`, ie. 48 index 0 has the function. 49 * `adj_inputs` is a list with the same number of items as 50 `inputs`, has the values of the adjoint input 51 * Other inputs unused. 52 53 Returns: the resulting product accounting for fixed DoFs 54 """ 55 adj_input = adj_inputs[0] 56 x = inputs[0].vector() 57 58 output = adj_input 59 # the derivative of the fixed dofs is null 60 output[self.fixed_indexes] = 0.0 61 return output
Zero-out the components of fixed DoFs for derivative.
inputs
is the list of inputs togel.fix_dofs.fix_dofs
, ie. index 0 has the function.adj_inputs
is a list with the same number of items asinputs
, has the values of the adjoint input- Other inputs unused.
Returns: the resulting product accounting for fixed DoFs
def
recompute_component(self, inputs, block_variable, idx, prepared):
63 def recompute_component(self, inputs, block_variable, idx, prepared): 64 """Re-fixes the DoFs 65 66 * `inputs` is the list of inputs to `fix_dofs` originally 67 encountered 68 * Rest of inputs unused 69 70 Returns: the DoF-fixed function. 71 """ 72 return _backend_fix_dofs( 73 inputs[0], 74 self.fixed_indexes, 75 self.fixed_values 76 )
Re-fixes the DoFs
inputs
is the list of inputs tofix_dofs
originally encountered- Rest of inputs unused
Returns: the DoF-fixed function.
def
fix_dofs(*args, **kwargs):
23 def _overloaded_function(*args, **kwargs): 24 annotate = annotate_tape(kwargs) 25 if annotate: 26 b_kwargs = block_class.pop_kwargs(kwargs) 27 b_kwargs.update(kwargs) 28 block = block_class(*args, **b_kwargs) 29 30 with stop_annotating(): 31 func_output = func(*args, **kwargs) 32 33 r = [] 34 func_output = Enlist(func_output) 35 for out in func_output: 36 r.append(create_overloaded_object(out)) 37 38 if annotate: 39 for out in r: 40 block.add_output(out.create_block_variable()) 41 tape = get_working_tape() 42 tape.add_block(block) 43 44 return func_output.delist(r)
Forces DoF values of a FE function to specific values.
Is equipped for automatic differentiation.
func
is a FEniCS function in a space with DoFsfixed_indexes
is an iterable of integer local DoF indicesfixed_values
is an iterable of the same size asfixed_indexes
with the float values to set the DoFs to.