Skip to content

Commit

Permalink
Add 3D case support capability for serial FSI cases (#133)
Browse files Browse the repository at this point in the history
* Initial un-tested design for 3D FEniCS cases

* autopep8 formatting checks

* Filter duplicate DoFs

* Filtering DoFs which are on the coupling interface

* Code refactoring to avoid parallel initialization for serial cases

* 3D function_space is identified as VECTOR FunctionType

* Rectify error in generating PointSource

* autopep8 formatting correction

* Remove unnecessary import

* Introducing changes which handle dimentionality in a generic way

* Add assertions to safeguard against use of FEniCS Expression type boundary conditions in 3D cases

* Checking if read and write functions passed are of the correct dimensions
  • Loading branch information
IshaanDesai authored Feb 7, 2022
1 parent 54c1642 commit a3f91d0
Show file tree
Hide file tree
Showing 2 changed files with 22 additions and 18 deletions.
16 changes: 10 additions & 6 deletions fenicsprecice/adapter_core.py
Original file line number Diff line number Diff line change
Expand Up @@ -78,24 +78,28 @@ class CouplingMode(Enum):
UNI_DIRECTIONAL_READ_COUPLING = 6


def determine_function_type(input_obj):
def determine_function_type(input_obj, dims):
"""
Determines if the function is scalar- or vector-valued based on rank evaluation.
Parameters
----------
input_obj :
A FEniCS function.
input_obj : FEniCS Function or FEniCS FunctionSpace
A FEniCS Function object or a FEniCS FunctionSpace object
dims : int
Dimension of problem.
Returns
-------
tag : bool
0 if input_function is SCALAR and 1 if input_function is VECTOR.
"""
if isinstance(input_obj, FunctionSpace): # scalar-valued functions have rank 0 is FEniCS
if input_obj.num_sub_spaces() == 0:
if isinstance(input_obj, FunctionSpace):
obj_dim = input_obj.num_sub_spaces()
if obj_dim == 0:
return FunctionType.SCALAR
elif input_obj.num_sub_spaces() == 2:
elif obj_dim == 2 or obj_dim == 3:
assert obj_dim == dims
return FunctionType.VECTOR
elif isinstance(input_obj, Function):
if input_obj.value_rank() == 0:
Expand Down
24 changes: 12 additions & 12 deletions fenicsprecice/fenicsprecice.py
Original file line number Diff line number Diff line change
Expand Up @@ -117,6 +117,7 @@ def create_coupling_expression(self):
coupling_expression : Object of class dolfin.functions.expression.Expression
Reference to object of class SegregatedRBFInterpolationExpression.
"""
assert(self._fenics_dims == 2), "Boundary conditions of Expression objects are only allowed for 2D cases"

if not (self._read_function_type is FunctionType.SCALAR or self._read_function_type is FunctionType.VECTOR):
raise Exception("No valid read_function is provided in initialization. Cannot create coupling expression")
Expand Down Expand Up @@ -159,6 +160,8 @@ def update_coupling_expression(self, coupling_expression, data):
The coupling data. A dictionary containing nodal data with vertex coordinates as key and associated data as
value.
"""
assert(self._fenics_dims == 2), "Boundary conditions of Expression objects are only allowed for 2D cases"

if not self._empty_rank:
coupling_expression.update_boundary_data(np.array(list(data.values())), np.array(list(data.keys())))

Expand Down Expand Up @@ -244,7 +247,7 @@ def write_data(self, write_function):
assert (w_func != write_function)

# Check that the function provided lives on the same function space provided during initialization
assert (self._write_function_type == determine_function_type(w_func))
assert (self._write_function_type == determine_function_type(w_func, self._fenics_dims))
assert (write_function.function_space() == self._write_function_space)

write_data_id = self._interface.get_data_id(self._config.get_write_data_name(),
Expand All @@ -254,7 +257,7 @@ def write_data(self, write_function):
assert (self._is_parallel()
), "having participants without coupling mesh nodes is only valid for parallel runs"

write_function_type = determine_function_type(write_function)
write_function_type = determine_function_type(write_function, self._fenics_dims)
assert (write_function_type in list(FunctionType))
write_data = convert_fenics_to_precice(write_function, self._owned_vertices.get_local_ids())
if write_function_type is FunctionType.SCALAR:
Expand Down Expand Up @@ -335,20 +338,20 @@ def initialize(self, coupling_subdomain, read_function_space=None, write_object=
raise Exception("Incorrect read and write function space combination provided. Please check input "
"parameters in initialization")

coords = function_space.tabulate_dof_coordinates()
_, self._fenics_dims = coords.shape

if self._coupling_type is CouplingMode.UNI_DIRECTIONAL_READ_COUPLING or \
self._coupling_type is CouplingMode.BI_DIRECTIONAL_COUPLING:
self._read_function_type = determine_function_type(read_function_space)
self._read_function_type = determine_function_type(read_function_space, self._fenics_dims)
self._read_function_space = read_function_space

if self._coupling_type is CouplingMode.UNI_DIRECTIONAL_WRITE_COUPLING or \
self._coupling_type is CouplingMode.BI_DIRECTIONAL_COUPLING:
# Ensure that function spaces of read and write functions are defined using the same mesh
self._write_function_type = determine_function_type(write_function_space)
self._write_function_type = determine_function_type(write_function_space, self._fenics_dims)
self._write_function_space = write_function_space

coords = function_space.tabulate_dof_coordinates()
_, self._fenics_dims = coords.shape

# Ensure that function spaces of read and write functions use the same mesh
if self._coupling_type is CouplingMode.BI_DIRECTIONAL_COUPLING:
assert (self._read_function_space.mesh() is write_function_space.mesh()
Expand All @@ -357,9 +360,6 @@ def initialize(self, coupling_subdomain, read_function_space=None, write_object=
if fixed_boundary:
self._Dirichlet_Boundary = fixed_boundary

if self._fenics_dims != 2:
raise Exception("Currently the fenics-adapter only supports 2D cases")

if self._fenics_dims != self._interface.get_dimensions():
raise Exception("Dimension of preCICE setup and FEniCS do not match")

Expand Down Expand Up @@ -393,8 +393,8 @@ def initialize(self, coupling_subdomain, read_function_space=None, write_object=
self._precice_vertex_ids = self._interface.set_mesh_vertices(self._interface.get_mesh_id(
self._config.get_coupling_mesh_name()), self._owned_vertices.get_coordinates())

if self._coupling_type is CouplingMode.UNI_DIRECTIONAL_READ_COUPLING or \
self._coupling_type is CouplingMode.BI_DIRECTIONAL_COUPLING:
if (self._coupling_type is CouplingMode.UNI_DIRECTIONAL_READ_COUPLING or
self._coupling_type is CouplingMode.BI_DIRECTIONAL_COUPLING) and self._is_parallel:
# Determine shared vertices with neighbouring processes and get dictionaries for communication
self._send_map, self._recv_map = get_communication_map(self._comm, self._read_function_space,
self._owned_vertices,
Expand Down

0 comments on commit a3f91d0

Please sign in to comment.