Coverage for  / dolfinx-env / lib / python3.12 / site-packages / io4dolfinx / checkpointing.py: 97%

207 statements  

« prev     ^ index     » next       coverage.py v7.14.0, created at 2026-05-12 11:21 +0000

1# Copyright (C) 2023-2026 Jørgen Schartum Dokken 

2# 

3# This file is part of io4dolfinx 

4# 

5# SPDX-License-Identifier: MIT 

6 

7from __future__ import annotations 

8 

9import typing 

10from pathlib import Path 

11from typing import Any 

12 

13from mpi4py import MPI 

14 

15import basix 

16import dolfinx 

17import numpy as np 

18import numpy.typing as npt 

19import ufl 

20from packaging.version import Version 

21 

22from . import compat 

23from .backends import FileMode, ReadMode, get_backend 

24from .comm_helpers import ( 

25 send_and_recv_cell_perm, 

26 send_dofmap_and_recv_values, 

27 send_dofs_and_recv_values, 

28) 

29from .readers import create_geometry_function_space 

30from .structures import ArrayData, FunctionData, MeshTagsData 

31from .utils import ( 

32 check_file_exists, 

33 compute_dofmap_pos, 

34 compute_local_range, 

35 index_owner, 

36 unroll_dofmap, 

37 unroll_insert_position, 

38) 

39from .writers import prepare_meshdata_for_storage 

40from .writers import write_function as _internal_function_writer 

41from .writers import write_mesh as _internal_mesh_writer 

42 

43__all__ = [ 

44 "read_mesh", 

45 "write_function", 

46 "read_function", 

47 "write_mesh", 

48 "read_meshtags", 

49 "write_meshtags", 

50 "read_attributes", 

51 "write_attributes", 

52] 

53 

54 

55def write_attributes( 

56 filename: Path | str, 

57 comm: MPI.Intracomm, 

58 name: str, 

59 attributes: dict[str, np.ndarray], 

60 backend_args: dict[str, typing.Any] | None = None, 

61 backend: str = "adios2", 

62): 

63 """Write attributes to file. 

64 

65 Args: 

66 filename: Path to file to write to 

67 comm: MPI communicator used in storage 

68 name: Name of the attributes 

69 attributes: Dictionary of attributes to write to file 

70 backend_args: Arguments for backend, for instance file type. 

71 backend: What backend to use for writing. 

72 """ 

73 backend_cls = get_backend(backend) 

74 backend_args = backend_cls.get_default_backend_args(backend_args) 

75 backend_cls.write_attributes(filename, comm, name, attributes, backend_args) 

76 

77 

78def read_attributes( 

79 filename: Path | str, 

80 comm: MPI.Intracomm, 

81 name: str, 

82 backend_args: dict[str, typing.Any] | None = None, 

83 backend: str = "adios2", 

84) -> dict[str, typing.Any]: 

85 """Read attributes from file. 

86 

87 Args: 

88 filename: Path to file to read from 

89 comm: MPI communicator used in storage 

90 name: Name of the attributes 

91 backend_args: Arguments for backend, for instance file type. 

92 backend: What backend to use for writing. 

93 Returns: 

94 The attributes 

95 """ 

96 backend_cls = get_backend(backend) 

97 backend_args = backend_cls.get_default_backend_args(backend_args) 

98 return backend_cls.read_attributes(filename, comm, name, backend_args) 

99 

100 

101def read_timestamps( 

102 filename: Path | str, 

103 comm: MPI.Intracomm, 

104 function_name: str, 

105 backend_args: dict[str, typing.Any] | None = None, 

106 backend: str = "adios2", 

107) -> npt.NDArray[np.float64 | str]: # type: ignore[type-var] 

108 """ 

109 Read time-stamps from a checkpoint file. 

110 

111 Args: 

112 comm: MPI communicator 

113 filename: Path to file 

114 function_name: Name of the function to read time-stamps for 

115 backend_args: Arguments for backend, for instance file type. 

116 backend: What backend to use for writing. 

117 Returns: 

118 The time-stamps 

119 """ 

120 check_file_exists(filename) 

121 backend_cls = get_backend(backend) 

122 backend_args = backend_cls.get_default_backend_args(backend_args) 

123 return backend_cls.read_timestamps(filename, comm, function_name, backend_args) 

124 

125 

126def write_meshtags( 

127 filename: Path | str, 

128 mesh: dolfinx.mesh.Mesh, 

129 meshtags: dolfinx.mesh.MeshTags, 

130 meshtag_name: typing.Optional[str] = None, 

131 backend_args: dict[str, Any] | None = None, 

132 backend: str = "adios2", 

133): 

134 """ 

135 Write meshtags associated with input mesh to file. 

136 

137 .. note:: 

138 For this checkpoint to work, the mesh must be written to file 

139 using :func:`write_mesh` before calling this function. 

140 

141 Args: 

142 filename: Path to save meshtags (with file-extension) 

143 mesh: The mesh associated with the meshtags 

144 meshtags: The meshtags to write to file 

145 meshtag_name: Name of the meshtag. If None, the meshtag name is used. 

146 backend_args: Option to IO backend. 

147 backend: IO backend 

148 """ 

149 

150 # Extract data from meshtags (convert to global geometry node indices for each entity) 

151 tag_entities = meshtags.indices 

152 dim = meshtags.dim 

153 num_tag_entities_local = mesh.topology.index_map(dim).size_local 

154 local_tag_entities = tag_entities[tag_entities < num_tag_entities_local] 

155 local_values = meshtags.values[: len(local_tag_entities)] 

156 

157 num_saved_tag_entities = len(local_tag_entities) 

158 local_start = mesh.comm.exscan(num_saved_tag_entities, op=MPI.SUM) 

159 local_start = local_start if mesh.comm.rank != 0 else 0 

160 global_num_tag_entities = mesh.comm.allreduce(num_saved_tag_entities, op=MPI.SUM) 

161 

162 dof_layout = compat.cmap(mesh).create_dof_layout() 

163 if hasattr(dof_layout, "num_entity_closure_dofs"): 

164 num_dofs_per_entity = dof_layout.num_entity_closure_dofs(dim) 

165 else: 

166 num_dofs_per_entity = len(dof_layout.entity_closure_dofs(dim, 0)) 

167 

168 entities_to_geometry = dolfinx.cpp.mesh.entities_to_geometry( 

169 mesh._cpp_object, dim, local_tag_entities, False 

170 ) 

171 

172 indices = ( 

173 mesh.geometry.index_map() 

174 .local_to_global(entities_to_geometry.reshape(-1)) 

175 .reshape(entities_to_geometry.shape) 

176 ) 

177 name = meshtag_name or meshtags.name 

178 

179 tag_ct = dolfinx.cpp.mesh.cell_entity_type(mesh.topology.cell_type, dim, 0).name 

180 tag_data = MeshTagsData( 

181 values=local_values, 

182 num_entities_global=global_num_tag_entities, 

183 num_dofs_per_entity=num_dofs_per_entity, 

184 indices=indices, 

185 name=name, 

186 local_start=local_start, 

187 dim=meshtags.dim, 

188 cell_type=tag_ct, 

189 ) 

190 

191 # Get backend and default arguments 

192 backend_cls = get_backend(backend) 

193 backend_args = backend_cls.get_default_backend_args(backend_args) 

194 return backend_cls.write_meshtags(filename, mesh.comm, tag_data, backend_args=backend_args) 

195 

196 

197def read_meshtags( 

198 filename: Path | str, 

199 mesh: dolfinx.mesh.Mesh, 

200 meshtag_name: str, 

201 backend_args: dict[str, Any] | None = None, 

202 backend: str = "adios2", 

203) -> dolfinx.mesh.MeshTags: 

204 """ 

205 Read meshtags from file and return a :class:`dolfinx.mesh.MeshTags` object. 

206 

207 Args: 

208 filename: Path to meshtags file (with file-extension) 

209 mesh: The mesh associated with the meshtags 

210 meshtag_name: The name of the meshtag to read 

211 engine: Adios2 Engine 

212 Returns: 

213 The meshtags 

214 """ 

215 check_file_exists(filename) 

216 backend_cls = get_backend(backend) 

217 backend_args = backend_cls.get_default_backend_args(backend_args) 

218 data = backend_cls.read_meshtags_data(filename, mesh.comm, meshtag_name, backend_args) 

219 

220 local_entities, local_values = dolfinx.io.distribute_entity_data( 

221 mesh, int(data.dim), data.indices, data.values 

222 ) 

223 mesh.topology.create_connectivity(data.dim, 0) 

224 mesh.topology.create_connectivity(data.dim, mesh.topology.dim) 

225 

226 adj = dolfinx.graph.adjacencylist(local_entities) 

227 

228 local_values = np.array(local_values, dtype=np.int32) 

229 

230 mt = dolfinx.mesh.meshtags_from_entities(mesh, int(data.dim), adj, local_values) 

231 mt.name = meshtag_name 

232 return mt 

233 

234 

235def read_function( 

236 filename: Path | str, 

237 u: dolfinx.fem.Function, 

238 time: float = 0.0, 

239 name: str | None = None, 

240 backend_args: dict[str, Any] | None = None, 

241 backend: str = "adios2", 

242): 

243 """ 

244 Read checkpoint from file and fill it into `u`. 

245 

246 Args: 

247 filename: Path to checkpoint 

248 u: Function to fill 

249 time: Time-stamp associated with checkpoint 

250 name: If not provided, `u.name` is used to search through the input file for the function 

251 """ 

252 check_file_exists(filename) 

253 

254 mesh = u.function_space.mesh 

255 comm = mesh.comm 

256 if name is None: 

257 name = u.name 

258 

259 # ----------------------Step 1--------------------------------- 

260 # Compute index of input cells and get cell permutation 

261 num_owned_cells = mesh.topology.index_map(mesh.topology.dim).size_local 

262 input_cells = mesh.topology.original_cell_index[:num_owned_cells] 

263 mesh.topology.create_entity_permutations() 

264 cell_perm = mesh.topology.get_cell_permutation_info()[:num_owned_cells] 

265 

266 # Compute mesh->input communicator 

267 # 1.1 Compute mesh->input communicator 

268 backend_cls = get_backend(backend) 

269 owners: npt.NDArray[np.int32] 

270 if backend_cls.read_mode == ReadMode.serial: 

271 owners = np.zeros(input_cells, dtype=np.int32) 

272 elif backend_cls.read_mode == ReadMode.parallel: 

273 num_cells_global = mesh.topology.index_map(mesh.topology.dim).size_global 

274 owners = index_owner(mesh.comm, input_cells, num_cells_global) 

275 else: 

276 raise NotImplementedError(f"{backend_cls.read_mode} not implemented") 

277 # -------------------Step 2------------------------------------ 

278 # Send and receive global cell index and cell perm 

279 inc_cells, inc_perms = send_and_recv_cell_perm(input_cells, cell_perm, owners, mesh.comm) 

280 

281 # -------------------Step 3----------------------------------- 

282 # Read dofmap from file and compute dof owners 

283 check_file_exists(filename) 

284 backend_cls = get_backend(backend) 

285 backend_args = backend_cls.get_default_backend_args(backend_args) 

286 

287 input_dofmap = backend_cls.read_dofmap(filename, comm, name, backend_args) 

288 

289 # Compute owner of dofs in dofmap 

290 dof_owner: npt.NDArray[np.int32] 

291 if backend_cls.read_mode == ReadMode.serial: 

292 dof_owner = np.zeros(len(input_dofmap.array), dtype=np.int32) 

293 elif backend_cls.read_mode == ReadMode.parallel: 

294 num_dofs_global = ( 

295 u.function_space.dofmap.index_map.size_global * u.function_space.dofmap.index_map_bs 

296 ) 

297 dof_owner = index_owner(comm, input_dofmap.array.astype(np.int64), num_dofs_global) 

298 else: 

299 raise NotImplementedError(f"{backend_cls.read_mode} not implemented") 

300 

301 # --------------------Step 4----------------------------------- 

302 # Read array from file and communicate them to input dofmap process 

303 input_array, starting_pos = backend_cls.read_dofs(filename, comm, name, time, backend_args) 

304 

305 recv_array = send_dofs_and_recv_values( 

306 input_dofmap.array.astype(np.int64), dof_owner, comm, input_array, starting_pos 

307 ) 

308 

309 # -------------------Step 5-------------------------------------- 

310 # Invert permutation of input data based on input perm 

311 # Then apply current permutation to the local data 

312 element = u.function_space.element 

313 if element.needs_dof_transformations: 

314 bs = u.function_space.dofmap.bs 

315 

316 # Read input cell permutations on dofmap process 

317 local_input_range = compute_local_range(comm, num_cells_global) 

318 input_local_cell_index = inc_cells - local_input_range[0] 

319 input_perms = backend_cls.read_cell_perms(comm, filename, backend_args) 

320 

321 # Start by sorting data array by cell permutation 

322 num_dofs_per_cell = input_dofmap.offsets[1:] - input_dofmap.offsets[:-1] 

323 assert np.allclose(num_dofs_per_cell, num_dofs_per_cell[0]) 

324 

325 # Sort dofmap by input local cell index 

326 input_perms_sorted = input_perms[input_local_cell_index] 

327 unrolled_dofmap_position = unroll_insert_position( 

328 input_local_cell_index, num_dofs_per_cell[0] 

329 ) 

330 dofmap_sorted_by_input = recv_array[unrolled_dofmap_position] 

331 

332 # First invert input data to reference element then transform to current mesh 

333 element.Tt_apply(dofmap_sorted_by_input, input_perms_sorted, bs) 

334 element.Tt_inv_apply(dofmap_sorted_by_input, inc_perms, bs) 

335 # Compute invert permutation 

336 inverted_perm = np.empty_like(unrolled_dofmap_position) 

337 inverted_perm[unrolled_dofmap_position] = np.arange( 

338 len(unrolled_dofmap_position), dtype=inverted_perm.dtype 

339 ) 

340 recv_array = dofmap_sorted_by_input[inverted_perm] 

341 

342 # ------------------Step 6---------------------------------------- 

343 # For each dof owned by a process, find the local position in the dofmap. 

344 V = u.function_space 

345 local_cells, dof_pos = compute_dofmap_pos(V) 

346 input_cells = V.mesh.topology.original_cell_index[local_cells] 

347 num_cells_global = V.mesh.topology.index_map(V.mesh.topology.dim).size_global 

348 

349 if backend_cls.read_mode == ReadMode.serial: 

350 owners = np.zeros(len(input_cells), dtype=np.int32) 

351 elif backend_cls.read_mode == ReadMode.parallel: 

352 owners = index_owner(V.mesh.comm, input_cells, num_cells_global) 

353 else: 

354 raise NotImplementedError(f"{backend_cls.read_mode} not implemented") 

355 

356 unique_owners, owner_count = np.unique(owners, return_counts=True) 

357 # FIXME: In C++ use NBX to find neighbourhood 

358 sub_comm = V.mesh.comm.Create_dist_graph( 

359 [V.mesh.comm.rank], [len(unique_owners)], unique_owners, reorder=False 

360 ) 

361 source, dest, _ = sub_comm.Get_dist_neighbors() 

362 sub_comm.Free() 

363 

364 owned_values = send_dofmap_and_recv_values( 

365 comm, 

366 np.asarray(source, dtype=np.int32), 

367 np.asarray(dest, dtype=np.int32), 

368 owners, 

369 owner_count.astype(np.int32), 

370 input_cells, 

371 dof_pos, 

372 num_cells_global, 

373 recv_array, 

374 input_dofmap.offsets, 

375 ) 

376 u.x.array[: len(owned_values)] = owned_values 

377 u.x.scatter_forward() 

378 

379 

380def read_mesh( 

381 filename: Path | str, 

382 comm: MPI.Intracomm, 

383 ghost_mode: dolfinx.mesh.GhostMode = dolfinx.mesh.GhostMode.shared_facet, 

384 time: float | str | None = 0.0, 

385 read_from_partition: bool = False, 

386 backend_args: dict[str, Any] | None = None, 

387 backend: str = "adios2", 

388 max_facet_to_cell_links: int = 2, 

389) -> dolfinx.mesh.Mesh: 

390 """ 

391 Read an ADIOS2 mesh into DOLFINx. 

392 

393 Args: 

394 filename: Path to input file 

395 comm: The MPI communciator to distribute the mesh over 

396 ghost_mode: Ghost mode to use for mesh. If `read_from_partition` 

397 is set to `True` this option is ignored. 

398 time: Time stamp associated with mesh 

399 read_from_partition: Read mesh with partition from file 

400 backend_args: List of arguments to reader backend 

401 max_facet_to_cell_links: Maximum number of cells a facet 

402 can be connected to. 

403 Returns: 

404 The distributed mesh 

405 """ 

406 # Read in data in a distributed fashin 

407 check_file_exists(filename) 

408 backend_cls = get_backend(backend) 

409 backend_args = backend_cls.get_default_backend_args(backend_args) 

410 

411 # Let each backend handle what should be default behavior when reading mesh 

412 # with or without time stamp. 

413 dist_in_data = backend_cls.read_mesh_data( 

414 filename, 

415 comm, 

416 time=time, 

417 read_from_partition=read_from_partition, 

418 backend_args=backend_args, 

419 ) 

420 

421 # Create DOLFINx mesh 

422 element = basix.ufl.element( 

423 basix.ElementFamily.P, 

424 dist_in_data.cell_type, 

425 dist_in_data.degree, 

426 basix.LagrangeVariant(int(dist_in_data.lvar)), 

427 shape=(dist_in_data.x.shape[1],), 

428 dtype=dist_in_data.x.dtype, 

429 ) 

430 domain = ufl.Mesh(element) 

431 if (partition_graph := dist_in_data.partition_graph) is not None: 

432 

433 def partitioner(comm: MPI.Intracomm, n, m, topo): 

434 assert len(topo[0]) % (len(partition_graph.offsets) - 1) == 0 

435 if Version(dolfinx.__version__) > Version("0.9.0"): 

436 return partition_graph._cpp_object 

437 else: 

438 return partition_graph 

439 else: 

440 try: 

441 partitioner = dolfinx.cpp.mesh.create_cell_partitioner( 

442 ghost_mode, max_facet_to_cell_links=max_facet_to_cell_links 

443 ) 

444 except TypeError: 

445 partitioner = dolfinx.cpp.mesh.create_cell_partitioner(ghost_mode) 

446 

447 # Should change to the commented code below when we require python 

448 # minimum version to be >=3.12 see https://github.com/python/cpython/pull/116198 

449 # import inspect 

450 # sig = inspect.signature(dolfinx.mesh.create_cell_partitioner) 

451 # part_kwargs = {} 

452 # if "max_facet_to_cell_links" in list(sig.parameters.keys()): 

453 # part_kwargs["max_facet_to_cell_links"] = max_facet_to_cell_links 

454 # partitioner = dolfinx.cpp.mesh.create_cell_partitioner(ghost_mode, **part_kwargs) 

455 

456 return dolfinx.mesh.create_mesh( 

457 comm, 

458 cells=dist_in_data.cells, 

459 x=dist_in_data.x, 

460 e=domain, 

461 partitioner=partitioner, 

462 ) 

463 

464 

465def write_mesh( 

466 filename: Path, 

467 mesh: dolfinx.mesh.Mesh, 

468 mode: FileMode = FileMode.write, 

469 time: float = 0.0, 

470 store_partition_info: bool = False, 

471 backend_args: dict[str, Any] | None = None, 

472 backend: str = "adios2", 

473): 

474 """ 

475 Write a mesh to file. 

476 

477 Args: 

478 filename: Path to save mesh (without file-extension) 

479 mesh: The mesh to write to file 

480 

481 store_partition_info: Store mesh partitioning (including ghosting) to file 

482 """ 

483 mesh_data = prepare_meshdata_for_storage(mesh=mesh, store_partition_info=store_partition_info) 

484 

485 _internal_mesh_writer( 

486 filename, 

487 mesh.comm, 

488 mesh_data=mesh_data, 

489 time=time, 

490 backend_args=backend_args, 

491 backend=backend, 

492 mode=mode, 

493 ) 

494 

495 

496def write_function( 

497 filename: Path | str, 

498 u: dolfinx.fem.Function, 

499 time: float = 0.0, 

500 mode: FileMode = FileMode.append, 

501 name: str | None = None, 

502 backend_args: dict[str, Any] | None = None, 

503 backend: str = "adios2", 

504): 

505 """ 

506 Write function checkpoint to file. 

507 

508 Args: 

509 u: Function to write to file 

510 time: Time-stamp for simulation 

511 filename: Path to write to 

512 mode: Write or append. 

513 name: Name of function to write. If None, the name of the function is used. 

514 backend_args: Arguments to the IO backend. 

515 backend: The backend to use 

516 """ 

517 dofmap = u.function_space.dofmap 

518 values = u.x.array 

519 mesh = u.function_space.mesh 

520 comm = mesh.comm 

521 mesh.topology.create_entity_permutations() 

522 cell_perm = mesh.topology.get_cell_permutation_info() 

523 num_cells_local = mesh.topology.index_map(mesh.topology.dim).size_local 

524 local_cell_range = mesh.topology.index_map(mesh.topology.dim).local_range 

525 num_cells_global = mesh.topology.index_map(mesh.topology.dim).size_global 

526 

527 # Convert local dofmap into global_dofmap 

528 dmap = dofmap.list 

529 num_dofs_per_cell = dmap.shape[1] 

530 dofmap_bs = dofmap.bs 

531 num_dofs_local_dmap = num_cells_local * num_dofs_per_cell * dofmap_bs 

532 index_map_bs = dofmap.index_map_bs 

533 

534 # Unroll dofmap for block size 

535 unrolled_dofmap = unroll_dofmap(dofmap.list[:num_cells_local, :], dofmap_bs) 

536 dmap_loc = (unrolled_dofmap // index_map_bs).reshape(-1) 

537 dmap_rem = (unrolled_dofmap % index_map_bs).reshape(-1) 

538 

539 # Convert imap index to global index 

540 imap_global = dofmap.index_map.local_to_global(dmap_loc) 

541 dofmap_global = imap_global * index_map_bs + dmap_rem 

542 dofmap_imap = dolfinx.common.IndexMap(mesh.comm, num_dofs_local_dmap) 

543 

544 # Compute dofmap offsets 

545 local_dofmap_offsets = np.arange(num_cells_local + 1, dtype=np.int64) 

546 local_dofmap_offsets[:] *= num_dofs_per_cell * dofmap_bs 

547 local_dofmap_offsets += dofmap_imap.local_range[0] 

548 

549 num_dofs_global = dofmap.index_map.size_global * dofmap.index_map_bs 

550 local_dof_range = np.asarray(dofmap.index_map.local_range) * dofmap.index_map_bs 

551 num_dofs_local = local_dof_range[1] - local_dof_range[0] 

552 

553 # Create internal data structure for function data to write to file 

554 function_data = FunctionData( 

555 cell_permutations=cell_perm[:num_cells_local].copy(), 

556 local_cell_range=local_cell_range, 

557 num_cells_global=num_cells_global, 

558 dofmap_array=dofmap_global, 

559 dofmap_offsets=local_dofmap_offsets, 

560 dofmap_range=dofmap_imap.local_range, 

561 global_dofs_in_dofmap=dofmap_imap.size_global, 

562 values=values[:num_dofs_local].copy(), 

563 dof_range=local_dof_range, 

564 num_dofs_global=num_dofs_global, 

565 name=name or u.name, 

566 ) 

567 # Write to file 

568 fname = Path(filename) 

569 _internal_function_writer( 

570 fname, comm, function_data, time, backend_args=backend_args, backend=backend, mode=mode 

571 ) 

572 

573 

574def read_function_names( 

575 filename: Path | str, 

576 comm: MPI.Intracomm, 

577 backend_args: dict[str, Any] | None = None, 

578 backend: str = "h5py", 

579) -> list[str]: 

580 """Read all function names from a file. 

581 

582 Args: 

583 filename: Path to file 

584 comm: MPI communicator to launch IO on. 

585 backend_args: Arguments to backend 

586 

587 Returns: 

588 A list of function names. 

589 """ 

590 backend_cls = get_backend(backend) 

591 return backend_cls.read_function_names(filename, comm, backend_args=backend_args) 

592 

593 

594def write_point_data( 

595 filename: Path | str, 

596 u: dolfinx.fem.Function, 

597 time: str | float | None, 

598 mode: FileMode, 

599 backend_args: dict[str, Any] | None, 

600 backend: str = "vtkhdf", 

601): 

602 """Write function to file by interpolating into geometry nodes. 

603 

604 

605 Args: 

606 filename: Path to file 

607 u: The function to store 

608 time: Time stamp 

609 mode: Append or write 

610 backend_args: The backend arguments 

611 backend: Which backend to use. 

612 """ 

613 V = create_geometry_function_space(u.function_space.mesh, int(np.prod(u.ufl_shape))) 

614 v_out = dolfinx.fem.Function(V, name=u.name, dtype=u.x.array.dtype) 

615 v_out.interpolate(u) 

616 comm = v_out.function_space.mesh.comm 

617 data_shape = (V.dofmap.index_map.size_global, V.dofmap.index_map_bs) 

618 local_range = V.dofmap.index_map.local_range 

619 num_dofs_local = V.dofmap.index_map.size_local 

620 data = v_out.x.array.reshape(-1, V.dofmap.index_map_bs)[:num_dofs_local] 

621 ad = ArrayData( 

622 name=v_out.name, values=data, global_shape=data_shape, local_range=local_range, type="Point" 

623 ) 

624 backend_cls = get_backend(backend) 

625 return backend_cls.write_data( 

626 filename, comm=comm, mode=mode, time=time, array_data=ad, backend_args=backend_args 

627 ) 

628 

629 

630def write_cell_data( 

631 filename: Path | str, 

632 u: dolfinx.fem.Function, 

633 time: str | float | None, 

634 mode: FileMode, 

635 backend_args: dict[str, Any] | None, 

636 backend: str = "vtkhdf", 

637): 

638 """Write function to file by interpolating into cell midpoints. 

639 

640 

641 Args: 

642 filename: Path to file 

643 point_data: Data to write to file 

644 time: Time stamp 

645 mode: Append or write 

646 backend_args: The backend arguments 

647 """ 

648 V = dolfinx.fem.functionspace(u.function_space.mesh, ("DG", 0, u.ufl_shape)) 

649 v_out = dolfinx.fem.Function(V, name=u.name, dtype=u.x.array.dtype) 

650 v_out.interpolate(u) 

651 comm = v_out.function_space.mesh.comm 

652 data_shape = (V.dofmap.index_map.size_global, V.dofmap.index_map_bs) 

653 local_range = V.dofmap.index_map.local_range 

654 num_dofs_local = V.dofmap.index_map.size_local 

655 data = v_out.x.array.reshape(-1, V.dofmap.index_map_bs)[:num_dofs_local] 

656 

657 backend_cls = get_backend(backend) 

658 ad = ArrayData( 

659 name=v_out.name, values=data, global_shape=data_shape, local_range=local_range, type="Cell" 

660 ) 

661 backend_cls = get_backend(backend) 

662 return backend_cls.write_data( 

663 filename, comm=comm, mode=mode, time=time, array_data=ad, backend_args=backend_args 

664 )