Neko 1.99.1
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-2023, 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!
35 use comm
36 use mpi_f08
40 use stl, only : stl_hdr_t, stl_triangle_t
41 implicit none
42 private
43
44 type(mpi_datatype) :: mpi_nmsh_hex
45 type(mpi_datatype) :: mpi_nmsh_quad
46 type(mpi_datatype) :: mpi_nmsh_zone
47 type(mpi_datatype) :: mpi_nmsh_curve
48
49 type(mpi_datatype) :: mpi_re2v1_data_xyz
50 type(mpi_datatype) :: mpi_re2v1_data_xy
51 type(mpi_datatype) :: mpi_re2v1_data_cv
52 type(mpi_datatype) :: mpi_re2v1_data_bc
53
54 type(mpi_datatype) :: mpi_re2v2_data_xyz
55 type(mpi_datatype) :: mpi_re2v2_data_xy
56 type(mpi_datatype) :: mpi_re2v2_data_cv
57 type(mpi_datatype) :: mpi_re2v2_data_bc
58
59 type(mpi_datatype) :: mpi_stl_header
60 type(mpi_datatype) :: mpi_stl_triangle
61
62 integer :: mpi_real_size
68
69 ! Public derived types and size definitions
80
81 ! Public subroutines
83
84contains
85
88 integer :: ierr
89
90 ! Define derived types
95
100
103
104 ! Check sizes of MPI types
105 call mpi_type_size(mpi_real, mpi_real_size, ierr)
106 call mpi_type_size(mpi_double_precision, mpi_double_precision_size, ierr)
107 call mpi_type_size(mpi_character, mpi_character_size, ierr)
108 call mpi_type_size(mpi_integer, mpi_integer_size, ierr)
109 call mpi_type_size(mpi_logical, mpi_logical_size, ierr)
110 call mpi_type_size(mpi_real_precision, mpi_real_prec_size, ierr)
111
112 call mpi_barrier(neko_comm, ierr)
113
114 end subroutine neko_mpi_types_init
115
118 type(nmsh_hex_t) :: nmsh_hex
119 type(mpi_datatype) :: type(17)
120 integer(kind=MPI_ADDRESS_KIND) :: disp(17), base
121 integer :: len(17), i, ierr
122
123 call mpi_get_address(nmsh_hex%el_idx, disp(1), ierr)
124 call mpi_get_address(nmsh_hex%v(1)%v_idx, disp(2), ierr)
125 call mpi_get_address(nmsh_hex%v(1)%v_xyz, disp(3), ierr)
126 call mpi_get_address(nmsh_hex%v(2)%v_idx, disp(4), ierr)
127 call mpi_get_address(nmsh_hex%v(2)%v_xyz, disp(5), ierr)
128 call mpi_get_address(nmsh_hex%v(3)%v_idx, disp(6), ierr)
129 call mpi_get_address(nmsh_hex%v(3)%v_xyz, disp(7), ierr)
130 call mpi_get_address(nmsh_hex%v(4)%v_idx, disp(8), ierr)
131 call mpi_get_address(nmsh_hex%v(4)%v_xyz, disp(9), ierr)
132 call mpi_get_address(nmsh_hex%v(5)%v_idx, disp(10), ierr)
133 call mpi_get_address(nmsh_hex%v(5)%v_xyz, disp(11), ierr)
134 call mpi_get_address(nmsh_hex%v(6)%v_idx, disp(12), ierr)
135 call mpi_get_address(nmsh_hex%v(6)%v_xyz, disp(13), ierr)
136 call mpi_get_address(nmsh_hex%v(7)%v_idx, disp(14), ierr)
137 call mpi_get_address(nmsh_hex%v(7)%v_xyz, disp(15), ierr)
138 call mpi_get_address(nmsh_hex%v(8)%v_idx, disp(16), ierr)
139 call mpi_get_address(nmsh_hex%v(8)%v_xyz, disp(17), ierr)
140
141
142 base = disp(1)
143 do i = 1, 17
144 disp(i) = mpi_aint_diff(disp(i), base)
145 end do
146
147 len(1) = 1
148 len(2:16:2) = 1
149 len(3:17:2) = 3
150
151 type(1) = mpi_integer
152 type(2:16:2) = mpi_integer
153 type(3:17:2) = mpi_double_precision
154 call mpi_type_create_struct(17, len, disp, type, mpi_nmsh_hex, ierr)
155 call mpi_type_commit(mpi_nmsh_hex, ierr)
156 end subroutine mpi_type_nmsh_hex_init
157
160 type(nmsh_quad_t) :: nmsh_quad
161 type(mpi_datatype) :: type(9)
162 integer(kind=MPI_ADDRESS_KIND) :: disp(9), base
163 integer :: len(9), i, ierr
164
165 call mpi_get_address(nmsh_quad%el_idx, disp(1), ierr)
166 call mpi_get_address(nmsh_quad%v(1)%v_idx, disp(2), ierr)
167 call mpi_get_address(nmsh_quad%v(1)%v_xyz, disp(3), ierr)
168 call mpi_get_address(nmsh_quad%v(2)%v_idx, disp(4), ierr)
169 call mpi_get_address(nmsh_quad%v(2)%v_xyz, disp(5), ierr)
170 call mpi_get_address(nmsh_quad%v(3)%v_idx, disp(6), ierr)
171 call mpi_get_address(nmsh_quad%v(3)%v_xyz, disp(7), ierr)
172 call mpi_get_address(nmsh_quad%v(4)%v_idx, disp(8), ierr)
173 call mpi_get_address(nmsh_quad%v(4)%v_xyz, disp(9), ierr)
174
175
176 base = disp(1)
177 do i = 1, 9
178 disp(i) = mpi_aint_diff(disp(i), base)
179 end do
180
181 len(1) = 1
182 len(2:8:2) = 1
183 len(3:9:2) = 3
184
185 type(1) = mpi_integer
186 type(2:8:2) = mpi_integer
187 type(3:9:2) = mpi_double_precision
188 call mpi_type_create_struct(9, len, disp, type, mpi_nmsh_quad, ierr)
189 call mpi_type_commit(mpi_nmsh_quad, ierr)
190 end subroutine mpi_type_nmsh_quad_init
191
194 type(nmsh_zone_t) :: nmsh_zone
195 type(mpi_datatype) :: type(6)
196 integer(kind=MPI_ADDRESS_KIND) :: disp(6), base
197 integer :: len(6), i, ierr
198
199 call mpi_get_address(nmsh_zone%e, disp(1), ierr)
200 call mpi_get_address(nmsh_zone%f, disp(2), ierr)
201 call mpi_get_address(nmsh_zone%p_e, disp(3), ierr)
202 call mpi_get_address(nmsh_zone%p_f, disp(4), ierr)
203 call mpi_get_address(nmsh_zone%glb_pt_ids, disp(5), ierr)
204 call mpi_get_address(nmsh_zone%type, disp(6), ierr)
205
206 base = disp(1)
207 do i = 1, 6
208 disp(i) = mpi_aint_diff(disp(i), base)
209 end do
210
211 len(1:4) = 1
212 len(5) = 4
213 len(6) = 1
214 type = mpi_integer
215
216 call mpi_type_create_struct(6, len, disp, type, mpi_nmsh_zone, ierr)
217 call mpi_type_commit(mpi_nmsh_zone, ierr)
218
219 end subroutine mpi_type_nmsh_zone_init
220
223 type(nmsh_curve_el_t) :: nmsh_curve_el
224 type(mpi_datatype) :: type(3)
225 integer(kind=MPI_ADDRESS_KIND) :: disp(3), base
226 integer :: len(3), i, ierr
227
228 call mpi_get_address(nmsh_curve_el%e, disp(1), ierr)
229 call mpi_get_address(nmsh_curve_el%curve_data, disp(2), ierr)
230 call mpi_get_address(nmsh_curve_el%type, disp(3), ierr)
231
232 base = disp(1)
233 do i = 1, 3
234 disp(i) = mpi_aint_diff(disp(i), base)
235 end do
236
237 len(1) = 1
238 len(2) = 5*12
239 len(3) = 12
240 type(1) = mpi_integer
241 type(2) = mpi_double_precision
242 type(3) = mpi_integer
243
244 call mpi_type_create_struct(3, len, disp, type, mpi_nmsh_curve, ierr)
245 call mpi_type_commit(mpi_nmsh_curve, ierr)
246
247 end subroutine mpi_type_nmsh_curve_init
248
249
252 type(re2v1_xyz_t) :: re2v1_data
253 type(re2v2_xyz_t) :: re2v2_data
254 type(mpi_datatype) :: type(4)
255 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
256 integer :: len(4), ierr, i
257
258 !
259 ! Setup version 1
260 !
261
262 call mpi_get_address(re2v1_data%rgroup, disp(1), ierr)
263 call mpi_get_address(re2v1_data%x, disp(2), ierr)
264 call mpi_get_address(re2v1_data%y, disp(3), ierr)
265 call mpi_get_address(re2v1_data%z, disp(4), ierr)
266
267 base = disp(1)
268 do i = 1, 4
269 disp(i) = mpi_aint_diff(disp(i), base)
270 end do
271
272 len = 8
273 len(1) = 1
274 type = mpi_real
275
276 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_xyz, ierr)
277 call mpi_type_commit(mpi_re2v1_data_xyz, ierr)
278
279 !
280 ! Setup version 2
281 !
282
283 call mpi_get_address(re2v2_data%rgroup, disp(1), ierr)
284 call mpi_get_address(re2v2_data%x, disp(2), ierr)
285 call mpi_get_address(re2v2_data%y, disp(3), ierr)
286 call mpi_get_address(re2v2_data%z, disp(4), ierr)
287
288 base = disp(1)
289 do i = 1, 4
290 disp(i) = mpi_aint_diff(disp(i), base)
291 end do
292
293 len = 8
294 len(1) = 1
295 type = mpi_double_precision
296
297 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_xyz, ierr)
298 call mpi_type_commit(mpi_re2v2_data_xyz, ierr)
299
300 end subroutine mpi_type_re2_xyz_init
301
304 type(re2v1_xy_t) :: re2v1_data
305 type(re2v2_xy_t) :: re2v2_data
306 type(mpi_datatype) :: type(3)
307 integer(kind=MPI_ADDRESS_KIND) :: disp(3), base
308 integer :: len(3), ierr, i
309
310 !
311 ! Setup version 1
312 !
313
314 call mpi_get_address(re2v1_data%rgroup, disp(1), ierr)
315 call mpi_get_address(re2v1_data%x, disp(2), ierr)
316 call mpi_get_address(re2v1_data%y, disp(3), ierr)
317
318 base = disp(1)
319 do i = 1, 3
320 disp(i) = mpi_aint_diff(disp(i), base)
321 end do
322
323 len = 4
324 len(1) = 1
325 type = mpi_real
326
327 call mpi_type_create_struct(3, len, disp, type, mpi_re2v1_data_xy, ierr)
328 call mpi_type_commit(mpi_re2v1_data_xy, ierr)
329
330 !
331 ! Setup version 2
332 !
333
334 call mpi_get_address(re2v2_data%rgroup, disp(1), ierr)
335 call mpi_get_address(re2v2_data%x, disp(2), ierr)
336 call mpi_get_address(re2v2_data%y, disp(3), ierr)
337
338 base = disp(1)
339 do i = 1, 3
340 disp(i) = mpi_aint_diff(disp(i), base)
341 end do
342
343 len = 4
344 len(1) = 1
345 type = mpi_double_precision
346
347 call mpi_type_create_struct(3, len, disp, type, mpi_re2v2_data_xy, ierr)
348 call mpi_type_commit(mpi_re2v2_data_xy, ierr)
349
350 end subroutine mpi_type_re2_xy_init
351
354 type(re2v1_curve_t) :: re2v1_data
355 type(re2v2_curve_t) :: re2v2_data
356 type(mpi_datatype) :: type(4)
357 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
358 integer :: len(4), ierr, i
359
360 !
361 ! Setup version 1
362 !
363
364 call mpi_get_address(re2v1_data%elem, disp(1), ierr)
365 call mpi_get_address(re2v1_data%zone, disp(2), ierr)
366 call mpi_get_address(re2v1_data%point, disp(3), ierr)
367 call mpi_get_address(re2v1_data%type, disp(4), ierr)
368
369 base = disp(1)
370 do i = 1, 4
371 disp(i) = mpi_aint_diff(disp(i), base)
372 end do
373
374 len(1:2) = 1
375 len(3) = 5
376 len(4) = 4
377 type(1:2) = mpi_integer
378 type(3) = mpi_real
379 type(4) = mpi_character
380
381 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_cv, ierr)
382 call mpi_type_commit(mpi_re2v1_data_cv, ierr)
383
384 !
385 ! Setup version 2
386 !
387
388 call mpi_get_address(re2v2_data%elem, disp(1), ierr)
389 call mpi_get_address(re2v2_data%zone, disp(2), ierr)
390 call mpi_get_address(re2v2_data%point, disp(3), ierr)
391 call mpi_get_address(re2v2_data%type, disp(4), ierr)
392
393 base = disp(1)
394 do i = 1, 4
395 disp(i) = mpi_aint_diff(disp(i), base)
396 end do
397
398 len(1:2) = 1
399 len(3) = 5
400 len(4) = 8
401 type(1:2) = mpi_double_precision
402 type(3) = mpi_double_precision
403 type(4) = mpi_character
404
405 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_cv, ierr)
406 call mpi_type_commit(mpi_re2v2_data_cv, ierr)
407
408 end subroutine mpi_type_re2_cv_init
409
412 type(re2v1_bc_t) :: re2v1_data
413 type(re2v2_bc_t) :: re2v2_data
414 type(mpi_datatype) :: type(4)
415 integer(kind=MPI_ADDRESS_KIND) :: disp(4), base
416 integer :: len(4), ierr, i
417
418 !
419 ! Setup version 1
420 !
421
422 call mpi_get_address(re2v1_data%elem, disp(1), ierr)
423 call mpi_get_address(re2v1_data%face, disp(2), ierr)
424 call mpi_get_address(re2v1_data%bc_data, disp(3), ierr)
425 call mpi_get_address(re2v1_data%type, disp(4), ierr)
426
427 base = disp(1)
428 do i = 1, 4
429 disp(i) = mpi_aint_diff(disp(i), base)
430 end do
431
432 len(1:2) = 1
433 len(3) = 5
434 len(4) = 4
435 type(1:2) = mpi_integer
436 type(3) = mpi_real
437 type(4) = mpi_character
438
439 call mpi_type_create_struct(4, len, disp, type, mpi_re2v1_data_bc, ierr)
440 call mpi_type_commit(mpi_re2v1_data_bc, ierr)
441
442 !
443 ! Setup version 2
444 !
445
446 call mpi_get_address(re2v2_data%elem, disp(1), ierr)
447 call mpi_get_address(re2v2_data%face, disp(2), ierr)
448 call mpi_get_address(re2v2_data%bc_data, disp(3), ierr)
449 call mpi_get_address(re2v2_data%type, disp(4), ierr)
450
451 base = disp(1)
452 do i = 1, 4
453 disp(i) = mpi_aint_diff(disp(i), base)
454 end do
455
456 len(1:2) = 1
457 len(3) = 5
458 len(4) = 8
459 type(1:2) = mpi_double_precision
460 type(3) = mpi_double_precision
461 type(4) = mpi_character
462
463 call mpi_type_create_struct(4, len, disp, type, mpi_re2v2_data_bc, ierr)
464 call mpi_type_commit(mpi_re2v2_data_bc, ierr)
465
466 end subroutine mpi_type_re2_bc_init
467
470 type(stl_hdr_t) :: stl_hdr
471 type(mpi_datatype) :: type(2)
472 integer(kind=MPI_ADDRESS_KIND) :: disp(2), base
473 integer :: len(2), ierr, i
474
475 call mpi_get_address(stl_hdr%hdr, disp(1), ierr)
476 call mpi_get_address(stl_hdr%ntri, disp(2), ierr)
477
478 base = disp(1)
479 do i = 1, 2
480 disp(i) = mpi_aint_diff(disp(i), base)
481 end do
482
483 len(1) = 80
484 len(2) = 1
485
486 type(1) = mpi_character
487 type(2) = mpi_integer
488
489 call mpi_type_create_struct(2, len, disp, type, mpi_stl_header, ierr)
490 call mpi_type_commit(mpi_stl_header, ierr)
491
492 end subroutine mpi_type_stl_header_init
493
496 type(stl_triangle_t) :: tri
497 type(mpi_datatype) :: type(5)
498 integer(kind=MPI_ADDRESS_KIND) :: disp(5), base
499 integer :: len(5), i, ierr
500
501 call mpi_get_address(tri%n, disp(1), ierr)
502 call mpi_get_address(tri%v1, disp(2), ierr)
503 call mpi_get_address(tri%v2, disp(3), ierr)
504 call mpi_get_address(tri%v3, disp(4), ierr)
505 call mpi_get_address(tri%attrib, disp(5), ierr)
506
507 base = disp(1)
508 do i = 1, 5
509 disp(i) = mpi_aint_diff(disp(i), base)
510 end do
511
512 len(1:4) = 3
513 len(5) = 1
514
515 type(1:4) = mpi_real
516 type(5) = mpi_integer2
517
518 call mpi_type_create_struct(5, len, disp, type, mpi_stl_triangle, ierr)
519 call mpi_type_commit(mpi_stl_triangle, ierr)
520
521 end subroutine mpi_type_stl_triangle_init
522
534 end subroutine neko_mpi_types_free
535
538 integer ierr
539 call mpi_type_free(mpi_nmsh_hex, ierr)
540 end subroutine mpi_type_nmsh_hex_free
541
544 integer ierr
545 call mpi_type_free(mpi_nmsh_quad, ierr)
546 end subroutine mpi_type_nmsh_quad_free
547
550 integer ierr
551 call mpi_type_free(mpi_nmsh_zone, ierr)
552 end subroutine mpi_type_nmsh_zone_free
553
556 integer ierr
557 call mpi_type_free(mpi_nmsh_curve, ierr)
558 end subroutine mpi_type_nmsh_curve_free
559
562 integer ierr
563 call mpi_type_free(mpi_re2v1_data_xyz, ierr)
564 call mpi_type_free(mpi_re2v2_data_xyz, ierr)
565 end subroutine mpi_type_re2_xyz_free
566
569 integer ierr
570 call mpi_type_free(mpi_re2v1_data_xy, ierr)
571 call mpi_type_free(mpi_re2v2_data_xy, ierr)
572 end subroutine mpi_type_re2_xy_free
573
576 integer ierr
577 call mpi_type_free(mpi_re2v1_data_cv, ierr)
578 call mpi_type_free(mpi_re2v2_data_cv, ierr)
579 end subroutine mpi_type_re2_cv_free
580
583 integer ierr
584 call mpi_type_free(mpi_re2v1_data_bc, ierr)
585 call mpi_type_free(mpi_re2v2_data_bc, ierr)
586 end subroutine mpi_type_re2_bc_free
587
590 integer ierr
591 call mpi_type_free(mpi_stl_header, ierr)
592 end subroutine mpi_type_stl_header_free
593
596 integer ierr
597 call mpi_type_free(mpi_stl_triangle, ierr)
598 end subroutine mpi_type_stl_triangle_free
599end 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:50
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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:44
subroutine, public neko_mpi_types_init
Define all MPI derived types.
Definition mpi_types.f90:88
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:63
integer, public mpi_real_prec_size
Size of working precision REAL types.
Definition mpi_types.f90:67
type(mpi_datatype), public mpi_re2v2_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:54
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:55
type(mpi_datatype), public mpi_nmsh_quad
MPI derived type for 2D Neko nmsh data.
Definition mpi_types.f90:45
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:47
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:64
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:62
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:65
type(mpi_datatype), public mpi_re2v1_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:51
type(mpi_datatype), public mpi_re2v2_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:56
integer, public mpi_logical_size
Size of MPI type logical.
Definition mpi_types.f90:66
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:50
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:52
type(mpi_datatype), public mpi_stl_header
MPI Derived type for a STL header.
Definition mpi_types.f90:59
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:46
type(mpi_datatype), public mpi_re2v2_data_bc
MPI derived type for NEKTON re2 bc data.
Definition mpi_types.f90:57
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:60
type(mpi_datatype), public mpi_re2v1_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:49
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