-
Notifications
You must be signed in to change notification settings - Fork 0
Read scalar values from OF #14
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Merged
Merged
Changes from all commits
Commits
Show all changes
6 commits
Select commit
Hold shift + click to select a range
807ce5b
added scalar mesh element and vector mesh element
RemDelaporteMathurin 526461c
zip file instead
RemDelaporteMathurin 0a7959f
Merge branch 'main' into read-scalar
RemDelaporteMathurin 40d3a6a
update with main
RemDelaporteMathurin 6dbc435
Merge branch 'main' into read-scalar
RemDelaporteMathurin dc5465c
dolfinx>=0.9
RemDelaporteMathurin File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -172,4 +172,6 @@ cython_debug/ | |
|
||
*.bp | ||
.vscode | ||
*_version.py | ||
*_version.py | ||
|
||
*Zone.Identifier |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -103,12 +103,15 @@ def _create_dolfinx_mesh(self): | |
) | ||
degree = 1 # Set polynomial degree | ||
cell = ufl.Cell(shape) | ||
self.mesh_element = basix.ufl.element( | ||
self.mesh_vector_element = basix.ufl.element( | ||
"Lagrange", cell.cellname(), degree, shape=(3,) | ||
) | ||
self.mesh_scalar_element = basix.ufl.element( | ||
"Lagrange", cell.cellname(), degree, shape=() | ||
) | ||
|
||
# Create dolfinx Mesh | ||
mesh_ufl = ufl.Mesh(self.mesh_element) | ||
mesh_ufl = ufl.Mesh(self.mesh_vector_element) | ||
self.dolfinx_mesh = create_mesh( | ||
MPI.COMM_WORLD, self.connectivity, self.OF_mesh.points, mesh_ufl | ||
) | ||
|
@@ -127,10 +130,17 @@ def create_dolfinx_function( | |
the dolfinx function | ||
""" | ||
self._read_with_pyvista(t=t) | ||
self._create_dolfinx_mesh() | ||
self.function_space = dolfinx.fem.functionspace( | ||
self.dolfinx_mesh, self.mesh_element | ||
) | ||
|
||
# Create dolfinx mesh if it doesn't exist | ||
if not hasattr(self, "dolfinx_mesh"): | ||
self._create_dolfinx_mesh() | ||
|
||
if name == "U": | ||
element = self.mesh_vector_element | ||
else: | ||
element = self.mesh_scalar_element | ||
Comment on lines
+138
to
+141
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Again could add small comment, in the case of reading for "U" a velocity field, a vector element is used etc |
||
|
||
self.function_space = dolfinx.fem.functionspace(self.dolfinx_mesh, element) | ||
u = dolfinx.fem.Function(self.function_space) | ||
|
||
num_vertices = ( | ||
|
Binary file not shown.
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -1,10 +1,35 @@ | ||
from foam2dolfinx import OpenFOAMReader | ||
import dolfinx | ||
from pyvista import examples | ||
import zipfile | ||
from pathlib import Path | ||
|
||
|
||
def test_reading_and_writing_cavity_example(): | ||
my_of_reader = OpenFOAMReader(filename=examples.download_cavity(load=False)) | ||
vel = my_of_reader.create_dolfinx_function(t=2.5) | ||
|
||
assert isinstance(vel, dolfinx.fem.Function) | ||
|
||
|
||
def test_baby_example(tmpdir): | ||
time = 2812.0 | ||
|
||
zip_path = Path("test/data/baby_example.zip") | ||
extract_path = Path(tmpdir) / "baby_example" | ||
|
||
# Unzip the file | ||
with zipfile.ZipFile(zip_path, "r") as zip_ref: | ||
zip_ref.extractall(extract_path) | ||
|
||
# Construct the path to the .foam file | ||
foam_file = extract_path / "baby_example/pv.foam" | ||
|
||
# read the .foam file | ||
my_of_reader = OpenFOAMReader(filename=str(foam_file), OF_mesh_cell_type_value=10) | ||
|
||
vel = my_of_reader.create_dolfinx_function(t=time, name="U") | ||
T = my_of_reader.create_dolfinx_function(t=time, name="T") | ||
|
||
assert isinstance(vel, dolfinx.fem.Function) | ||
assert isinstance(T, dolfinx.fem.Function) |
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Could probably add a small inline comment here stating how the type of element for creating the mesh is not important