Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
mpi_types.f90
Go to the documentation of this file.
1! Copyright (c) 2019-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
36 use mpi_f08, only : mpi_type_size, mpi_type_create_struct, mpi_type_commit, &
37 mpi_get_address, mpi_real, mpi_double_precision, mpi_character, &
38 mpi_integer, mpi_logical, mpi_integer2, mpi_address_kind, &
39 mpi_datatype, mpi_aint_diff
43 use stl, only : stl_hdr_t, stl_triangle_t
44 implicit none
45 private
46
47 type(mpi_datatype) :: mpi_nmsh_hex
48 type(mpi_datatype) :: mpi_nmsh_quad
49 type(mpi_datatype) :: mpi_nmsh_zone
50 type(mpi_datatype) :: mpi_nmsh_curve
51
52 type(mpi_datatype) :: mpi_re2v1_data_xyz
53 type(mpi_datatype) :: mpi_re2v1_data_xy
54 type(mpi_datatype) :: mpi_re2v1_data_cv
55 type(mpi_datatype) :: mpi_re2v1_data_bc
56
57 type(mpi_datatype) :: mpi_re2v2_data_xyz
58 type(mpi_datatype) :: mpi_re2v2_data_xy
59 type(mpi_datatype) :: mpi_re2v2_data_cv
60 type(mpi_datatype) :: mpi_re2v2_data_bc
61
62 type(mpi_datatype) :: mpi_stl_header
63 type(mpi_datatype) :: mpi_stl_triangle
64
65 integer :: mpi_real_size
71
72 ! Public derived types and size definitions
83
84 ! Public subroutines
86
87contains
88
91 integer :: ierr
92
93 ! Define derived types
98
103
106
107 ! Check sizes of MPI types
108 call mpi_type_size(mpi_real, mpi_real_size, ierr)
109 call mpi_type_size(mpi_double_precision, mpi_double_precision_size, ierr)
110 call mpi_type_size(mpi_character, mpi_character_size, ierr)
111 call mpi_type_size(mpi_integer, mpi_integer_size, ierr)
112 call mpi_type_size(mpi_logical, mpi_logical_size, ierr)
113 call mpi_type_size(mpi_real_precision, mpi_real_prec_size, ierr)
114
115 call mpi_barrier(neko_comm, ierr)
116
117 end subroutine neko_mpi_types_init
118
121 type(nmsh_hex_t) :: nmsh_hex
122 type(mpi_datatype) :: type(17)
123 integer(kind=MPI_ADDRESS_KIND) :: disp(17), base
124 integer :: len(17), i, ierr
125
126 call mpi_get_address(nmsh_hex%el_idx, disp(1), ierr)
127 call mpi_get_address(nmsh_hex%v(1)%v_idx, disp(2), ierr)
128 call mpi_get_address(nmsh_hex%v(1)%v_xyz, disp(3), ierr)
129 call mpi_get_address(nmsh_hex%v(2)%v_idx, disp(4), ierr)
130 call mpi_get_address(nmsh_hex%v(2)%v_xyz, disp(5), ierr)
131 call mpi_get_address(nmsh_hex%v(3)%v_idx, disp(6), ierr)
132 call mpi_get_address(nmsh_hex%v(3)%v_xyz, disp(7), ierr)
133 call mpi_get_address(nmsh_hex%v(4)%v_idx, disp(8), ierr)
134 call mpi_get_address(nmsh_hex%v(4)%v_xyz, disp(9), ierr)
135 call mpi_get_address(nmsh_hex%v(5)%v_idx, disp(10), ierr)
136 call mpi_get_address(nmsh_hex%v(5)%v_xyz, disp(11), ierr)
137 call mpi_get_address(nmsh_hex%v(6)%v_idx, disp(12), ierr)
138 call mpi_get_address(nmsh_hex%v(6)%v_xyz, disp(13), ierr)
139 call mpi_get_address(nmsh_hex%v(7)%v_idx, disp(14), ierr)
140 call mpi_get_address(nmsh_hex%v(7)%v_xyz, disp(15), ierr)
141 call mpi_get_address(nmsh_hex%v(8)%v_idx, disp(16), ierr)
142 call mpi_get_address(nmsh_hex%v(8)%v_xyz, disp(17), ierr)
143
144
145 base = disp(1)
146 do i = 1, 17
147 disp(i) = mpi_aint_diff(disp(i), base)
148 end do
149
150 len(1) = 1
151 len(2:16:2) = 1
152 len(3:17:2) = 3
153
154 type(1) = mpi_integer
155 type(2:16:2) = mpi_integer
156 type(3:17:2) = mpi_double_precision
157 call mpi_type_create_struct(17, len, disp, type, mpi_nmsh_hex, ierr)
158 call mpi_type_commit(mpi_nmsh_hex, ierr)
159 end subroutine mpi_type_nmsh_hex_init
160
163 type(nmsh_quad_t) :: nmsh_quad
164 type(mpi_datatype) :: type(9)
165 integer(kind=MPI_ADDRESS_KIND) :: disp(9), base
166 integer :: len(9), i, ierr
167
168 call mpi_get_address(nmsh_quad%el_idx, disp(1), ierr)
169 call mpi_get_address(nmsh_quad%v(1)%v_idx, disp(2), ierr)
170 call mpi_get_address(nmsh_quad%v(1)%v_xyz, disp(3), ierr)
171 call mpi_get_address(nmsh_quad%v(2)%v_idx, disp(4), ierr)
172 call mpi_get_address(nmsh_quad%v(2)%v_xyz, disp(5), ierr)
173 call mpi_get_address(nmsh_quad%v(3)%v_idx, disp(6), ierr)
174 call mpi_get_address(nmsh_quad%v(3)%v_xyz, disp(7), ierr)
175 call mpi_get_address(nmsh_quad%v(4)%v_idx, disp(8), ierr)
176 call mpi_get_address(nmsh_quad%v(4)%v_xyz, disp(9), ierr)
177
178
179 base = disp(1)
180 do i = 1, 9
181 disp(i) = mpi_aint_diff(disp(i), base)
182 end do
183
184 len(1) = 1
185 len(2:8:2) = 1
186 len(3:9:2) = 3
187
188 type(1) = mpi_integer
189 type(2:8:2) = mpi_integer
190 type(3:9:2) = mpi_double_precision
191 call mpi_type_create_struct(9, len, disp, type, mpi_nmsh_quad, ierr)
192 call mpi_type_commit(mpi_nmsh_quad, ierr)
193 end subroutine mpi_type_nmsh_quad_init
194
197 type(nmsh_zone_t) :: nmsh_zone
198 type(mpi_datatype) :: type(6)
199 integer(kind=MPI_ADDRESS_KIND) :: disp(6), base
200 integer :: len(6), i, ierr
201
202 call mpi_get_address(nmsh_zone%e, disp(1), ierr)
203 call mpi_get_address(nmsh_zone%f, disp(2), ierr)
204 call mpi_get_address(nmsh_zone%p_e, disp(3), ierr)
205 call mpi_get_address(nmsh_zone%p_f, disp(4), ierr)
206 call mpi_get_address(nmsh_zone%glb_pt_ids, disp(5), ierr)
207 call mpi_get_address(nmsh_zone%type, disp(6), ierr)
208
209 base = disp(1)
210 do i = 1, 6
211 disp(i) = mpi_aint_diff(disp(i), base)
212 end do
213
214 len(1:4) = 1
215 len(5) = 4
216 len(6) = 1
217 type = mpi_integer
218
219 call mpi_type_create_struct(6, len, disp, type, mpi_nmsh_zone, ierr)
220 call mpi_type_commit(mpi_nmsh_zone, ierr)
221
222 end subroutine mpi_type_nmsh_zone_init
223
226 type(nmsh_curve_el_t) :: nmsh_curve_el
227 type(mpi_datatype) :: type(3)
228 integer(kind=MPI_ADDRESS_KIND) :: disp(3), base
229 integer :: len(3), i, ierr
230
231 call mpi_get_address(nmsh_curve_el%e, disp(1), ierr)
232 call mpi_get_address(nmsh_curve_el%curve_data, disp(2), ierr)
233 call mpi_get_address(nmsh_curve_el%type, disp(3), ierr)
234
235 base = disp(1)
236 do i = 1, 3
237 disp(i) = mpi_aint_diff(disp(i), base)
238 end do
239
240 len(1) = 1
241 len(2) = 5*12
242 len(3) = 12
243 type(1) = mpi_integer
244 type(2) = mpi_double_precision
245 type(3) = mpi_integer
246
247 call mpi_type_create_struct(3, len, disp, type, mpi_nmsh_curve, ierr)
248 call mpi_type_commit(mpi_nmsh_curve, ierr)
249
250 end subroutine mpi_type_nmsh_curve_init
251
252
255 type(re2v1_xyz_t) :: re2v1_data
256 type(re2v2_xyz_t) :: re2v2_data
257 type(mpi_datatype) :: type(4)
258 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
259 integer :: len(4), ierr, i
260
261 !
262 ! Setup version 1
263 !
264
265 call mpi_get_address(re2v1_data%rgroup, disp(1), ierr)
266 call mpi_get_address(re2v1_data%x, disp(2), ierr)
267 call mpi_get_address(re2v1_data%y, disp(3), ierr)
268 call mpi_get_address(re2v1_data%z, disp(4), ierr)
269
270 base = disp(1)
271 do i = 1, 4
272 disp(i) = mpi_aint_diff(disp(i), base)
273 end do
274
275 len = 8
276 len(1) = 1
277 type = mpi_real
278
279 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_xyz, ierr)
280 call mpi_type_commit(mpi_re2v1_data_xyz, ierr)
281
282 !
283 ! Setup version 2
284 !
285
286 call mpi_get_address(re2v2_data%rgroup, disp(1), ierr)
287 call mpi_get_address(re2v2_data%x, disp(2), ierr)
288 call mpi_get_address(re2v2_data%y, disp(3), ierr)
289 call mpi_get_address(re2v2_data%z, disp(4), ierr)
290
291 base = disp(1)
292 do i = 1, 4
293 disp(i) = mpi_aint_diff(disp(i), base)
294 end do
295
296 len = 8
297 len(1) = 1
298 type = mpi_double_precision
299
300 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_xyz, ierr)
301 call mpi_type_commit(mpi_re2v2_data_xyz, ierr)
302
303 end subroutine mpi_type_re2_xyz_init
304
307 type(re2v1_xy_t) :: re2v1_data
308 type(re2v2_xy_t) :: re2v2_data
309 type(mpi_datatype) :: type(3)
310 integer(kind=MPI_ADDRESS_KIND) :: disp(3), base
311 integer :: len(3), ierr, i
312
313 !
314 ! Setup version 1
315 !
316
317 call mpi_get_address(re2v1_data%rgroup, disp(1), ierr)
318 call mpi_get_address(re2v1_data%x, disp(2), ierr)
319 call mpi_get_address(re2v1_data%y, disp(3), ierr)
320
321 base = disp(1)
322 do i = 1, 3
323 disp(i) = mpi_aint_diff(disp(i), base)
324 end do
325
326 len = 4
327 len(1) = 1
328 type = mpi_real
329
330 call mpi_type_create_struct(3, len, disp, type, mpi_re2v1_data_xy, ierr)
331 call mpi_type_commit(mpi_re2v1_data_xy, ierr)
332
333 !
334 ! Setup version 2
335 !
336
337 call mpi_get_address(re2v2_data%rgroup, disp(1), ierr)
338 call mpi_get_address(re2v2_data%x, disp(2), ierr)
339 call mpi_get_address(re2v2_data%y, disp(3), ierr)
340
341 base = disp(1)
342 do i = 1, 3
343 disp(i) = mpi_aint_diff(disp(i), base)
344 end do
345
346 len = 4
347 len(1) = 1
348 type = mpi_double_precision
349
350 call mpi_type_create_struct(3, len, disp, type, mpi_re2v2_data_xy, ierr)
351 call mpi_type_commit(mpi_re2v2_data_xy, ierr)
352
353 end subroutine mpi_type_re2_xy_init
354
357 type(re2v1_curve_t) :: re2v1_data
358 type(re2v2_curve_t) :: re2v2_data
359 type(mpi_datatype) :: type(4)
360 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
361 integer :: len(4), ierr, i
362
363 !
364 ! Setup version 1
365 !
366
367 call mpi_get_address(re2v1_data%elem, disp(1), ierr)
368 call mpi_get_address(re2v1_data%zone, disp(2), ierr)
369 call mpi_get_address(re2v1_data%point, disp(3), ierr)
370 call mpi_get_address(re2v1_data%type, disp(4), ierr)
371
372 base = disp(1)
373 do i = 1, 4
374 disp(i) = mpi_aint_diff(disp(i), base)
375 end do
376
377 len(1:2) = 1
378 len(3) = 5
379 len(4) = 4
380 type(1:2) = mpi_integer
381 type(3) = mpi_real
382 type(4) = mpi_character
383
384 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_cv, ierr)
385 call mpi_type_commit(mpi_re2v1_data_cv, ierr)
386
387 !
388 ! Setup version 2
389 !
390
391 call mpi_get_address(re2v2_data%elem, disp(1), ierr)
392 call mpi_get_address(re2v2_data%zone, disp(2), ierr)
393 call mpi_get_address(re2v2_data%point, disp(3), ierr)
394 call mpi_get_address(re2v2_data%type, disp(4), ierr)
395
396 base = disp(1)
397 do i = 1, 4
398 disp(i) = mpi_aint_diff(disp(i), base)
399 end do
400
401 len(1:2) = 1
402 len(3) = 5
403 len(4) = 8
404 type(1:2) = mpi_double_precision
405 type(3) = mpi_double_precision
406 type(4) = mpi_character
407
408 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_cv, ierr)
409 call mpi_type_commit(mpi_re2v2_data_cv, ierr)
410
411 end subroutine mpi_type_re2_cv_init
412
415 type(re2v1_bc_t) :: re2v1_data
416 type(re2v2_bc_t) :: re2v2_data
417 type(mpi_datatype) :: type(4)
418 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
419 integer :: len(4), ierr, i
420
421 !
422 ! Setup version 1
423 !
424
425 call mpi_get_address(re2v1_data%elem, disp(1), ierr)
426 call mpi_get_address(re2v1_data%face, disp(2), ierr)
427 call mpi_get_address(re2v1_data%bc_data, disp(3), ierr)
428 call mpi_get_address(re2v1_data%type, disp(4), ierr)
429
430 base = disp(1)
431 do i = 1, 4
432 disp(i) = mpi_aint_diff(disp(i), base)
433 end do
434
435 len(1:2) = 1
436 len(3) = 5
437 len(4) = 4
438 type(1:2) = mpi_integer
439 type(3) = mpi_real
440 type(4) = mpi_character
441
442 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_bc, ierr)
443 call mpi_type_commit(mpi_re2v1_data_bc, ierr)
444
445 !
446 ! Setup version 2
447 !
448
449 call mpi_get_address(re2v2_data%elem, disp(1), ierr)
450 call mpi_get_address(re2v2_data%face, disp(2), ierr)
451 call mpi_get_address(re2v2_data%bc_data, disp(3), ierr)
452 call mpi_get_address(re2v2_data%type, disp(4), ierr)
453
454 base = disp(1)
455 do i = 1, 4
456 disp(i) = mpi_aint_diff(disp(i), base)
457 end do
458
459 len(1:2) = 1
460 len(3) = 5
461 len(4) = 8
462 type(1:2) = mpi_double_precision
463 type(3) = mpi_double_precision
464 type(4) = mpi_character
465
466 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_bc, ierr)
467 call mpi_type_commit(mpi_re2v2_data_bc, ierr)
468
469 end subroutine mpi_type_re2_bc_init
470
473 type(stl_hdr_t) :: stl_hdr
474 type(mpi_datatype) :: type(2)
475 integer(kind=MPI_ADDRESS_KIND) :: disp(2), base
476 integer :: len(2), ierr, i
477
478 call mpi_get_address(stl_hdr%hdr, disp(1), ierr)
479 call mpi_get_address(stl_hdr%ntri, disp(2), ierr)
480
481 base = disp(1)
482 do i = 1, 2
483 disp(i) = mpi_aint_diff(disp(i), base)
484 end do
485
486 len(1) = 80
487 len(2) = 1
488
489 type(1) = mpi_character
490 type(2) = mpi_integer
491
492 call mpi_type_create_struct(2, len, disp, type, mpi_stl_header, ierr)
493 call mpi_type_commit(mpi_stl_header, ierr)
494
495 end subroutine mpi_type_stl_header_init
496
499 type(stl_triangle_t) :: tri
500 type(mpi_datatype) :: type(5)
501 integer(kind=MPI_ADDRESS_KIND) :: disp(5), base
502 integer :: len(5), i, ierr
503
504 call mpi_get_address(tri%n, disp(1), ierr)
505 call mpi_get_address(tri%v1, disp(2), ierr)
506 call mpi_get_address(tri%v2, disp(3), ierr)
507 call mpi_get_address(tri%v3, disp(4), ierr)
508 call mpi_get_address(tri%attrib, disp(5), ierr)
509
510 base = disp(1)
511 do i = 1, 5
512 disp(i) = mpi_aint_diff(disp(i), base)
513 end do
514
515 len(1:4) = 3
516 len(5) = 1
517
518 type(1:4) = mpi_real
519 type(5) = mpi_integer2
520
521 call mpi_type_create_struct(5, len, disp, type, mpi_stl_triangle, ierr)
522 call mpi_type_commit(mpi_stl_triangle, ierr)
523
524 end subroutine mpi_type_stl_triangle_init
525
537 end subroutine neko_mpi_types_free
538
541 integer ierr
542 call mpi_type_free(mpi_nmsh_hex, ierr)
543 end subroutine mpi_type_nmsh_hex_free
544
547 integer ierr
548 call mpi_type_free(mpi_nmsh_quad, ierr)
549 end subroutine mpi_type_nmsh_quad_free
550
553 integer ierr
554 call mpi_type_free(mpi_nmsh_zone, ierr)
555 end subroutine mpi_type_nmsh_zone_free
556
559 integer ierr
560 call mpi_type_free(mpi_nmsh_curve, ierr)
561 end subroutine mpi_type_nmsh_curve_free
562
565 integer ierr
566 call mpi_type_free(mpi_re2v1_data_xyz, ierr)
567 call mpi_type_free(mpi_re2v2_data_xyz, ierr)
568 end subroutine mpi_type_re2_xyz_free
569
572 integer ierr
573 call mpi_type_free(mpi_re2v1_data_xy, ierr)
574 call mpi_type_free(mpi_re2v2_data_xy, ierr)
575 end subroutine mpi_type_re2_xy_free
576
579 integer ierr
580 call mpi_type_free(mpi_re2v1_data_cv, ierr)
581 call mpi_type_free(mpi_re2v2_data_cv, ierr)
582 end subroutine mpi_type_re2_cv_free
583
586 integer ierr
587 call mpi_type_free(mpi_re2v1_data_bc, ierr)
588 call mpi_type_free(mpi_re2v2_data_bc, ierr)
589 end subroutine mpi_type_re2_bc_free
590
593 integer ierr
594 call mpi_type_free(mpi_stl_header, ierr)
595 end subroutine mpi_type_stl_header_free
596
599 integer ierr
600 call mpi_type_free(mpi_stl_triangle, ierr)
601 end subroutine mpi_type_stl_triangle_free
602end module neko_mpi_types
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
MPI derived types.
Definition mpi_types.f90:34
subroutine mpi_type_nmsh_quad_free
Deallocate nmsh quad derived MPI type.
subroutine mpi_type_re2_bc_free
Deallocate re2 bc derived MPI type.
subroutine mpi_type_nmsh_quad_init
Define a MPI derived type for a 2d nmsh quad.
type(mpi_datatype), public mpi_nmsh_hex
MPI derived type for 3D Neko nmsh data.
Definition mpi_types.f90:47
subroutine, public neko_mpi_types_init
Define all MPI derived types.
Definition mpi_types.f90:91
subroutine mpi_type_stl_header_init
Define a MPI dervied type for a STL header.
subroutine mpi_type_nmsh_curve_init
Define a MPI derived type for a nmsh curved element.
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition mpi_types.f90:66
integer, public mpi_real_prec_size
Size of working precision REAL types.
Definition mpi_types.f90:70
type(mpi_datatype), public mpi_re2v2_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:57
subroutine mpi_type_nmsh_zone_init
Define a MPI derived type for a nmsh zone.
type(mpi_datatype), public mpi_re2v2_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition mpi_types.f90:58
type(mpi_datatype), public mpi_nmsh_quad
MPI derived type for 2D Neko nmsh data.
Definition mpi_types.f90:48
subroutine mpi_type_stl_triangle_init
Define a MPI derived type for a STL triangle.
subroutine mpi_type_nmsh_hex_init
Define a MPI derived type for a 3d nmsh hex.
type(mpi_datatype), public mpi_nmsh_curve
MPI derived type for Neko nmsh curved elements.
Definition mpi_types.f90:50
subroutine mpi_type_re2_xyz_free
Deallocate re2 xyz derived MPI type.
subroutine mpi_type_stl_triangle_free
Deallocate STL triangle derived MPI type.
integer, public mpi_character_size
Size of MPI type character.
Definition mpi_types.f90:67
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:65
subroutine mpi_type_re2_cv_init
Define a MPI derived type for re2 cv data.
subroutine mpi_type_re2_xy_free
Deallocate re2 xyz derived MPI type.
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:68
type(mpi_datatype), public mpi_re2v1_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:54
type(mpi_datatype), public mpi_re2v2_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:59
integer, public mpi_logical_size
Size of MPI type logical.
Definition mpi_types.f90:69
subroutine mpi_type_nmsh_zone_free
Deallocate nmsh zone derived MPI type.
type(mpi_datatype), public mpi_re2v1_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition mpi_types.f90:53
subroutine mpi_type_nmsh_curve_free
Deallocate nmsh curve derived MPI type.
subroutine mpi_type_re2_bc_init
Define a MPI derived type for re2 bc data.
subroutine mpi_type_re2_cv_free
Deallocate re2 cv derived MPI type.
type(mpi_datatype), public mpi_re2v1_data_bc
MPI derived type for NEKTON re2 bc data.
Definition mpi_types.f90:55
type(mpi_datatype), public mpi_stl_header
MPI Derived type for a STL header.
Definition mpi_types.f90:62
subroutine mpi_type_nmsh_hex_free
Deallocate nmsh hex derived MPI type.
type(mpi_datatype), public mpi_nmsh_zone
MPI derived type for Neko nmsh zone data.
Definition mpi_types.f90:49
type(mpi_datatype), public mpi_re2v2_data_bc
MPI derived type for NEKTON re2 bc data.
Definition mpi_types.f90:60
subroutine mpi_type_re2_xy_init
Define a MPI derived type for a 2d re2 data.
subroutine mpi_type_re2_xyz_init
Define a MPI derived type for a 3d re2 data.
type(mpi_datatype), public mpi_stl_triangle
MPI derived type for a STL triangle.
Definition mpi_types.f90:63
type(mpi_datatype), public mpi_re2v1_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:52
subroutine, public neko_mpi_types_free
Deallocate all derived MPI types.
subroutine mpi_type_stl_header_free
Deallocate STL header dervied MPI type.
Neko binary mesh format.
Definition nmsh.f90:2
NEKTON re2 format.
Definition re2.f90:2
Stereolithography format.
Definition stl.f90:2
Defines a triangular element.
Definition tri.f90:34
Neko curve data.
Definition nmsh.f90:39
Neko hex element data.
Definition nmsh.f90:24
Neko quad element data.
Definition nmsh.f90:19
Neko zone data.
Definition nmsh.f90:29
NEKTON re2 bc data (version 1)
Definition re2.f90:39
NEKTON re2 curve data (version 1)
Definition re2.f90:31
NEKTON re2 element data (2d) (version 1)
Definition re2.f90:25
NEKTON re2 element data (3d) (version 1)
Definition re2.f90:18
NEKTON re2 bc data (version 2)
Definition re2.f90:73
NEKTON re2 curve data (version 2)
Definition re2.f90:65
NEKTON re2 element data (2d) (version 2)
Definition re2.f90:59
NEKTON re2 element data (3d) (version 2)
Definition re2.f90:52
Defines a STL hdr.
Definition stl.f90:8
Defines a STL triangle.
Definition stl.f90:14