ComIn 0.5.1
ICON Community Interface
Loading...
Searching...
No Matches
comin_descrdata.F90
Go to the documentation of this file.
1!> @file comin_descrdata.F90
2!! @brief Accessor functions for ComIn descriptive data structures.
3!
4! @authors 10/2021 :: ICON Community Interface <comin@icon-model.org>
5!
6! SPDX-License-Identifier: BSD-3-Clause
7!
8! See LICENSES for license information.
9! Where software is supplied by third parties, it is indicated in the
10! headers of the routines.
11!
13
14 use, INTRINSIC :: iso_c_binding, only: c_int, c_ptr, c_loc, c_null_char, c_loc, c_signed_char
15 USE comin_setup_constants, ONLY: wp
16 USE comin_state, ONLY: state
23
24 IMPLICIT NONE
25
26 PRIVATE
27
43
44 PUBLIC :: comin_descrdata_to_c
45
46#include "comin_global.inc"
47
48 INTERFACE is_in_bounds
49 MODULE PROCEDURE is_in_bounds_int
50 MODULE PROCEDURE is_in_bounds_ddd
51 END INTERFACE is_in_bounds
52
54 MODULE PROCEDURE comin_descrdata_to_c_signed_char_1
55 MODULE PROCEDURE comin_descrdata_to_c_string
56 MODULE PROCEDURE comin_descrdata_to_c_int_1
57 MODULE PROCEDURE comin_descrdata_to_c_int_2
58 MODULE PROCEDURE comin_descrdata_to_c_int_3
59 MODULE PROCEDURE comin_descrdata_to_c_real_1
60 MODULE PROCEDURE comin_descrdata_to_c_real_2
61 MODULE PROCEDURE comin_descrdata_to_c_real_3
62 MODULE PROCEDURE comin_descrdata_to_c_real_4
63 END INTERFACE comin_descrdata_to_c
64
65 INTERFACE
66 SUBROUTINE comin_current_set_datetime_c(datetime_c) BIND(C, NAME="comin_current_set_datetime")
67 USE iso_c_binding, ONLY: c_char
68 CHARACTER(kind=c_char), DIMENSION(*), INTENT(IN) :: datetime_c
69 END SUBROUTINE comin_current_set_datetime_c
70
71 FUNCTION comin_current_get_datetime_c() &
72 & result(datetime) &
73 & BIND(C, NAME="comin_current_get_datetime")
74 USE iso_c_binding, ONLY: c_ptr
75 TYPE(C_PTR) :: datetime
76 END FUNCTION comin_current_get_datetime_c
77
78 !> Conversion of global cell index to MPI-process local index.
79 !! @ingroup fortran_interface
80 INTEGER(C_INT) FUNCTION comin_descrdata_index_lookup_glb2loc_cell(jg, global_idx) &
81 & result(loc) BIND(C)
82 IMPORT :: c_int
83 INTEGER(kind=C_INT), INTENT(IN), VALUE :: jg !< domain index
84 INTEGER(kind=C_INT), INTENT(IN), VALUE :: global_idx !< global cell index
86
87 !> Conversion of global edge index to MPI-process local index.
88 !! @ingroup fortran_interface
89 INTEGER(C_INT) FUNCTION comin_descrdata_index_lookup_glb2loc_edge(jg, global_idx) &
90 & result(loc) BIND(C)
91 IMPORT :: c_int
92 INTEGER(kind=C_INT), INTENT(IN), VALUE :: jg !< domain index
93 INTEGER(kind=C_INT), INTENT(IN), VALUE :: global_idx !< global edge index
95
96 !> Conversion of global vert index to MPI-process local index.
97 !! @ingroup fortran_interface
98 INTEGER(C_INT) FUNCTION comin_descrdata_index_lookup_glb2loc_vert(jg, global_idx) &
99 & result(loc) BIND(C)
100 IMPORT :: c_int
101 INTEGER(kind=C_INT), INTENT(IN), VALUE :: jg !< domain index
102 INTEGER(kind=C_INT), INTENT(IN), VALUE :: global_idx !< global vert index
104
106 & RESULT(exp_start) &
107 & BIND(C)
108 use, INTRINSIC :: iso_c_binding, only: c_ptr
110
112 & RESULT(exp_stop) &
113 & BIND(C)
114 use, INTRINSIC :: iso_c_binding, only: c_ptr
116
118 & RESULT(run_start) &
119 & BIND(C)
120 use, INTRINSIC :: iso_c_binding, only: c_ptr
122
124 & RESULT(run_stop) &
125 & BIND(C)
126 use, INTRINSIC :: iso_c_binding, only: c_ptr
128
129 !> Receive pointer on array storing timestep information for all domains
130 !! @ingroup fortran_interface
132 use, INTRINSIC :: iso_c_binding, only: c_double, c_int
133 REAL(c_double) :: comin_descrdata_get_timesteplength
134 INTEGER(C_INT), INTENT(IN), VALUE :: jg
136
137 !> Fill array with timestep.
138 SUBROUTINE comin_descrdata_set_timesteplength(jg, dt_current) BIND(C)
139 use, INTRINSIC :: iso_c_binding, only: c_int, c_double
140 INTEGER(C_INT), INTENT(IN), VALUE :: jg
141 REAL(C_DOUBLE), INTENT(IN), VALUE :: dt_current
143
144 SUBROUTINE comin_state_set_simulation_interval(exp_start, exp_stop, run_start, run_stop) &
145 BIND(C)
146 use, INTRINSIC :: iso_c_binding, only: c_char
147 CHARACTER(KIND=C_CHAR, LEN=1), DIMENSION(*), INTENT(IN) :: exp_start, exp_stop, run_start, run_stop
149
150 END INTERFACE
151
152CONTAINS
153
154 !> Fill global data.
155 SUBROUTINE comin_descrdata_set_global(comin_global_info)
156 TYPE(t_comin_descrdata_global_data), INTENT(IN) :: comin_global_info
157 state%comin_descrdata_global_data = comin_global_info
158 END SUBROUTINE comin_descrdata_set_global
159
160 !> Set up data type for grid data.
161 SUBROUTINE comin_descrdata_set_domain(comin_domain_info)
162 TYPE(t_comin_descrdata_domain_data), INTENT(IN) :: comin_domain_info(:)
163 ALLOCATE(state%comin_descrdata_domain_data, source=comin_domain_info)
164 END SUBROUTINE comin_descrdata_set_domain
165
166 !> Fill time stamp info.
167 SUBROUTINE comin_descrdata_set_simulation_interval(exp_start, exp_stop, run_start, run_stop)
168 CHARACTER(LEN=*), INTENT(IN) :: exp_start, exp_stop, run_start, run_stop
169
171 trim(exp_start) // c_null_char, &
172 trim(exp_stop) // c_null_char, &
173 trim(run_start) // c_null_char, &
174 trim(run_stop) // c_null_char)
176
177 !> Retrieve time stamp info, current time information
178 !! @ingroup fortran_interface
179 SUBROUTINE comin_current_get_datetime(sim_time_current)
180 CHARACTER(LEN=:), ALLOCATABLE, INTENT(OUT) :: sim_time_current
181 !
182 TYPE(c_ptr) :: datetime_c
183 datetime_c = comin_current_get_datetime_c()
184 sim_time_current = convert_c_string(datetime_c)
185 END SUBROUTINE comin_current_get_datetime
186
187 !> Update time stamp info, current time information.
188 SUBROUTINE comin_current_set_datetime(sim_time_current)
189 CHARACTER(LEN=*), INTENT(IN) :: sim_time_current
190 CALL comin_current_set_datetime_c(trim(sim_time_current)//c_null_char)
191 END SUBROUTINE comin_current_set_datetime
192
193 !> request a pointer to simulation status
194 !! @ingroup fortran_interface
207
208 !> Clean descriptive data structure in ComIn
209 !> currently no content but keep for future use
211
212 END SUBROUTINE comin_descrdata_finalize
213
214 !!
215 !> auxiliary functions taken from ICON, version 2.6.5
216 !!
217
218 !> from mo_parallel_config
219 !> names of routines in ICON are: blk_no, idx_no
220 !-------------------------------------------------------------------------
221 ! The following two functions are for conversion of 1D to 2D indices and vice versa
222 !
223 ! Treatment of 0 (important for empty domains) and negative numbers:
224 !
225 ! Converting 1D => 2D:
226 !
227 ! 0 always is mapped to blk_no = 1, idx_no = 0
228 ! negative numbers: Convert usings ABS(j) and negate idx_no
229 !
230 ! Thus: blk_no >= 1 always!
231 ! idx_no > 0 for j > 0
232 ! idx_no = 0 for j = 0
233 ! idx_no < 0 for j < 0
234 !
235 ! This mimics mostly the behaviour of reshape_idx in mo_model_domimp_patches
236 ! with a difference for nproma=1 and j=0 (where reshape_idx returns blk_no=0, idx_no=1)
237 !
238 ! The consistent treatment of 0 in the above way is very important for empty domains
239 ! where start_index=1, end_index=0
240 !
241 ! Converting 2D => 1D:
242 ! Trying to invert the above and catching cases with blk_no < 1
243 !-------------------------------------------------------------------------
244
245 !> Auxiliary function: conversion of 1D to 2D indices.
246 !! @ingroup fortran_interface
247 INTEGER(c_int) FUNCTION comin_descrdata_get_block(idx1D) BIND(C, name="comin_descrdata_get_block")
248 INTEGER(c_int), INTENT(IN), VALUE :: idx1d
249 ! i.e. also 1 for idx1D=0, nproma=1
250 comin_descrdata_get_block = max((abs(idx1d)-1)/state%comin_descrdata_global_data%nproma + 1, 1)
251 END FUNCTION comin_descrdata_get_block
252
253 !> Auxiliary function: conversion of 1D to 2D indices.
254 !! @ingroup fortran_interface
255 INTEGER(c_int) FUNCTION comin_descrdata_get_index(idx1D) BIND(C, name="comin_descrdata_get_index")
256 INTEGER(c_int), INTENT(IN), VALUE :: idx1d
257 IF(idx1d==0) THEN
259 ELSE
260 comin_descrdata_get_index = sign(mod(abs(idx1d)-1,state%comin_descrdata_global_data%nproma)+1, idx1d)
261 ENDIF
262 END FUNCTION comin_descrdata_get_index
263
264 !> Computes the start and end block indices for loops over cell-based variables.
265 !! @ingroup fortran_interface
266 !!
267 SUBROUTINE comin_descrdata_get_cell_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk) &
268 & BIND(C, NAME="comin_descrdata_get_cell_block_limits")
269
270 INTEGER(c_int), INTENT(IN), VALUE :: jg !< Patch index for comin_domain.
271 INTEGER(c_int), INTENT(IN), VALUE :: irl_start !< `refin_ctrl` start block level.
272 INTEGER(c_int), INTENT(IN), VALUE :: irl_end !< `refin_ctrl` end block level.
273
274 INTEGER(c_int), INTENT(OUT) :: i_startblk !< Start index (jb loop)
275 INTEGER(c_int), INTENT(OUT) :: i_endblk !< End index (jb loop)
276
277 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
278 IF ( &
279 & is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%cells%start_block) .AND. &
280 & is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%cells%end_block)) THEN
281
282 i_startblk = state%comin_descrdata_domain_data(jg)%cells%start_block(irl_start)
283 i_endblk = state%comin_descrdata_domain_data(jg)%cells%end_block(irl_end)
284
285 ELSE
286 i_startblk = 1
287 i_endblk = 0
288 CALL comin_error_set(comin_error_index_out_of_bounds)
289 END IF
290 ELSE
291 i_startblk = 1
292 i_endblk = 0
293 CALL comin_error_set(comin_error_index_out_of_bounds)
294 END IF
295
297
298 !> Computes the start and end indices of do loops for cell-based variables.
299 !! @ingroup fortran_interface
300 !!
301 !! From ICON's `mo_loopindices`; name of corresponding ICON routine: `get_indices_c`.
302 !!
303 SUBROUTINE comin_descrdata_get_cell_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, &
304 i_endidx, irl_start, irl_end) &
305 & BIND(C, NAME="comin_descrdata_get_cell_indices")
306
307 INTEGER(c_int), INTENT(IN), VALUE :: jg ! Patch index for comin_domain
308 INTEGER(c_int), INTENT(IN), VALUE :: i_blk ! Current block (variable jb in do loops)
309 INTEGER(c_int), INTENT(IN), VALUE :: i_startblk ! Start block of do loop
310 INTEGER(c_int), INTENT(IN), VALUE :: i_endblk ! End block of do loop
311 INTEGER(c_int), INTENT(IN), VALUE :: irl_start ! refin_ctrl level where do loop starts
312 INTEGER(c_int), INTENT(IN), VALUE :: irl_end ! refin_ctrl level where do loop ends
313
314 INTEGER(c_int), INTENT(OUT) :: i_startidx, i_endidx ! Start and end indices (jc loop)
315
316 IF (i_blk == i_startblk) THEN
317 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
318 IF (is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%cells%start_block)) THEN
319 i_startidx = max(1,state%comin_descrdata_domain_data(jg)%cells%start_index(irl_start))
320 i_endidx = state%comin_descrdata_global_data%nproma
321 ELSE
322 i_startidx = 1
323 i_endidx = 0
324 CALL comin_error_set(comin_error_index_out_of_bounds)
325 END IF
326
327 IF (i_blk == i_endblk) THEN
328 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%cells%end_block)) THEN
329 i_endidx = state%comin_descrdata_domain_data(jg)%cells%end_index(irl_end)
330 ELSE
331 i_startidx = 1
332 i_endidx = 0
333 CALL comin_error_set(comin_error_index_out_of_bounds)
334 END IF
335 END IF
336 ELSE
337 i_startidx = 1
338 i_endidx = 0
339 CALL comin_error_set(comin_error_index_out_of_bounds)
340 END IF
341 ELSE IF (i_blk == i_endblk) THEN
342 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
343 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%cells%end_block)) THEN
344 i_startidx = 1
345 i_endidx = state%comin_descrdata_domain_data(jg)%cells%end_index(irl_end)
346 ELSE
347 i_startidx = 1
348 i_endidx = 0
349 CALL comin_error_set(comin_error_index_out_of_bounds)
350 END IF
351 ELSE
352 i_startidx = 1
353 i_endidx = 0
354 CALL comin_error_set(comin_error_index_out_of_bounds)
355 END IF
356 ELSE
357 i_startidx = 1
358 i_endidx = state%comin_descrdata_global_data%nproma
359 ENDIF
360
362
363 !> Computes the start and end block indices for loops over edge-based variables.
364 !! @ingroup fortran_interface
365 !!
366 SUBROUTINE comin_descrdata_get_edge_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk) &
367 & BIND(C, NAME="comin_descrdata_get_edge_block_limits")
368
369 INTEGER(c_int), INTENT(IN), VALUE :: jg !< Patch index for comin_domain.
370 INTEGER(c_int), INTENT(IN), VALUE :: irl_start !< `refin_ctrl` start block level.
371 INTEGER(c_int), INTENT(IN), VALUE :: irl_end !< `refin_ctrl` end block level.
372
373 INTEGER(c_int), INTENT(OUT) :: i_startblk !< Start index (jb loop)
374 INTEGER(c_int), INTENT(OUT) :: i_endblk !< End index (jb loop)
375
376 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
377 IF ( &
378 & is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%edges%start_block) .AND. &
379 & is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%edges%end_block)) THEN
380
381 i_startblk = state%comin_descrdata_domain_data(jg)%edges%start_block(irl_start)
382 i_endblk = state%comin_descrdata_domain_data(jg)%edges%end_block(irl_end)
383
384 ELSE
385 i_startblk = 1
386 i_endblk = 0
387 CALL comin_error_set(comin_error_index_out_of_bounds)
388 END IF
389 ELSE
390 i_startblk = 1
391 i_endblk = 0
392 CALL comin_error_set(comin_error_index_out_of_bounds)
393 END IF
394
396
397 !> Computes the start and end indices of do loops for edge-based variables.
398 !! @ingroup fortran_interface
399 !!
400 !! From ICON's `mo_loopindices`; name of corresponding ICON routine: `get_indices_e`.
401 !!
402 SUBROUTINE comin_descrdata_get_edge_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, &
403 i_endidx, irl_start, irl_end) &
404 & BIND(C, NAME="comin_descrdata_get_edge_indices")
405
406 INTEGER(c_int), INTENT(IN), VALUE :: jg ! Patch index for comin_domain
407 INTEGER(c_int), INTENT(IN), VALUE :: i_blk ! Current block (variable jb in do loops)
408 INTEGER(c_int), INTENT(IN), VALUE :: i_startblk ! Start block of do loop
409 INTEGER(c_int), INTENT(IN), VALUE :: i_endblk ! End block of do loop
410 INTEGER(c_int), INTENT(IN), VALUE :: irl_start ! refin_ctrl level where do loop starts
411 INTEGER(c_int), INTENT(IN), VALUE :: irl_end ! refin_ctrl level where do loop ends
412
413 INTEGER(c_int), INTENT(OUT) :: i_startidx, i_endidx ! Start and end indices (jc loop)
414
415 IF (i_blk == i_startblk) THEN
416 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
417 IF (is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%edges%start_block)) THEN
418 i_startidx = max(1,state%comin_descrdata_domain_data(jg)%edges%start_index(irl_start))
419 i_endidx = state%comin_descrdata_global_data%nproma
420 ELSE
421 i_startidx = 1
422 i_endidx = 0
423 CALL comin_error_set(comin_error_index_out_of_bounds)
424 END IF
425
426 IF (i_blk == i_endblk) THEN
427 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%edges%end_block)) THEN
428 i_endidx = state%comin_descrdata_domain_data(jg)%edges%end_index(irl_end)
429 ELSE
430 i_startidx = 1
431 i_endidx = 0
432 CALL comin_error_set(comin_error_index_out_of_bounds)
433 END IF
434 END IF
435 ELSE
436 i_startidx = 1
437 i_endidx = 0
438 CALL comin_error_set(comin_error_index_out_of_bounds)
439 END IF
440 ELSE IF (i_blk == i_endblk) THEN
441 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
442 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%edges%end_block)) THEN
443 i_startidx = 1
444 i_endidx = state%comin_descrdata_domain_data(jg)%edges%end_index(irl_end)
445 ELSE
446 i_startidx = 1
447 i_endidx = 0
448 CALL comin_error_set(comin_error_index_out_of_bounds)
449 END IF
450 ELSE
451 i_startidx = 1
452 i_endidx = 0
453 CALL comin_error_set(comin_error_index_out_of_bounds)
454 END IF
455 ELSE
456 i_startidx = 1
457 i_endidx = state%comin_descrdata_global_data%nproma
458 ENDIF
459
461
462 !> Computes the start and end block indices for loops over vertex-based variables.
463 !! @ingroup fortran_interface
464 !!
465 SUBROUTINE comin_descrdata_get_vert_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk) &
466 & BIND(C, NAME="comin_descrdata_get_vert_block_limits")
467
468 INTEGER(c_int), INTENT(IN), VALUE :: jg !< Patch index for comin_domain.
469 INTEGER(c_int), INTENT(IN), VALUE :: irl_start !< `refin_ctrl` start block level.
470 INTEGER(c_int), INTENT(IN), VALUE :: irl_end !< `refin_ctrl` end block level.
471
472 INTEGER(c_int), INTENT(OUT) :: i_startblk !< Start index (jb loop)
473 INTEGER(c_int), INTENT(OUT) :: i_endblk !< End index (jb loop)
474
475 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
476 IF ( &
477 & is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%verts%start_block) .AND. &
478 & is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%verts%end_block)) THEN
479
480 i_startblk = state%comin_descrdata_domain_data(jg)%verts%start_block(irl_start)
481 i_endblk = state%comin_descrdata_domain_data(jg)%verts%end_block(irl_end)
482
483 ELSE
484 i_startblk = 1
485 i_endblk = 0
486 CALL comin_error_set(comin_error_index_out_of_bounds)
487 END IF
488 ELSE
489 i_startblk = 1
490 i_endblk = 0
491 CALL comin_error_set(comin_error_index_out_of_bounds)
492 END IF
493
495
496 !> Computes the start and end indices of do loops for vertex-based variables.
497 !! @ingroup fortran_interface
498 !!
499 !! From ICON's `mo_loopindices`; name of corresponding ICON routine: `get_indices_v`.
500 !!
501 SUBROUTINE comin_descrdata_get_vert_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, &
502 i_endidx, irl_start, irl_end) &
503 & BIND(C, NAME="comin_descrdata_get_vert_indices")
504
505 INTEGER(c_int), INTENT(IN), VALUE :: jg ! Patch index for comin_domain
506 INTEGER(c_int), INTENT(IN), VALUE :: i_blk ! Current block (variable jb in do loops)
507 INTEGER(c_int), INTENT(IN), VALUE :: i_startblk ! Start block of do loop
508 INTEGER(c_int), INTENT(IN), VALUE :: i_endblk ! End block of do loop
509 INTEGER(c_int), INTENT(IN), VALUE :: irl_start ! refin_ctrl level where do loop starts
510 INTEGER(c_int), INTENT(IN), VALUE :: irl_end ! refin_ctrl level where do loop ends
511
512 INTEGER(c_int), INTENT(OUT) :: i_startidx, i_endidx ! Start and end indices (jc loop)
513
514 IF (i_blk == i_startblk) THEN
515 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
516 IF (is_in_bounds(irl_start, state%comin_descrdata_domain_data(jg)%verts%start_block)) THEN
517 i_startidx = max(1,state%comin_descrdata_domain_data(jg)%verts%start_index(irl_start))
518 i_endidx = state%comin_descrdata_global_data%nproma
519 ELSE
520 i_startidx = 1
521 i_endidx = 0
522 CALL comin_error_set(comin_error_index_out_of_bounds)
523 END IF
524
525 IF (i_blk == i_endblk) THEN
526 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%verts%end_block)) THEN
527 i_endidx = state%comin_descrdata_domain_data(jg)%verts%end_index(irl_end)
528 ELSE
529 i_startidx = 1
530 i_endidx = 0
531 CALL comin_error_set(comin_error_index_out_of_bounds)
532 END IF
533 END IF
534 ELSE
535 i_startidx = 1
536 i_endidx = 0
537 CALL comin_error_set(comin_error_index_out_of_bounds)
538 END IF
539 ELSE IF (i_blk == i_endblk) THEN
540 IF (is_in_bounds(jg, state%comin_descrdata_domain_data)) THEN
541 IF (is_in_bounds(irl_end, state%comin_descrdata_domain_data(jg)%verts%end_block)) THEN
542 i_startidx = 1
543 i_endidx = state%comin_descrdata_domain_data(jg)%verts%end_index(irl_end)
544 ELSE
545 i_startidx = 1
546 i_endidx = 0
547 CALL comin_error_set(comin_error_index_out_of_bounds)
548 END IF
549 ELSE
550 i_startidx = 1
551 i_endidx = 0
552 CALL comin_error_set(comin_error_index_out_of_bounds)
553 END IF
554 ELSE
555 i_startidx = 1
556 i_endidx = state%comin_descrdata_global_data%nproma
557 ENDIF
558
560
561 !> Calculate `npromz` value for the blocking, needed for patch allocation.
562 !> ... for the cells
563 !! @ingroup fortran_interface
564 !!
565 !! NB: Avoid the case nblks=0 for empty patches, this might cause troubles
566 !! if a empty patch is used somewhere (and npromz gets wrong in the formulas below).
567 !!
568 !! from `mo_setup_subdivision`; name of ICON routine: `npromz_c`.
569 INTEGER(c_int) FUNCTION comin_descrdata_get_cell_npromz(jg) BIND(C)
570 INTEGER(c_int), INTENT(IN), VALUE :: jg ! domain index for comin_domain
571
572 comin_descrdata_get_cell_npromz = state%comin_descrdata_domain_data(jg)%cells%ncells - &
573 & (state%comin_descrdata_domain_data(jg)%cells%nblks-1)*state%comin_descrdata_global_data%nproma
575
576 !> Calculate `npromz` value for the blocking, needed for patch allocation.
577 !> ... for the edges
578 !! @ingroup fortran_interface
579 !!
580 !! NB: Avoid the case nblks=0 for empty patches, this might cause troubles
581 !! if a empty patch is used somewhere (and npromz gets wrong in the formulas below).
582 !!
583 !! from `mo_setup_subdivision`; name of ICON routine: `npromz_e`.
584 INTEGER(c_int) FUNCTION comin_descrdata_get_edge_npromz(jg) BIND(C)
585 INTEGER(c_int), INTENT(IN), VALUE :: jg ! domain index for comin_domain
586
587 comin_descrdata_get_edge_npromz = state%comin_descrdata_domain_data(jg)%edges%nedges - &
588 & (state%comin_descrdata_domain_data(jg)%edges%nblks-1)*state%comin_descrdata_global_data%nproma
590
591 !> Calculate `npromz` value for the blocking, needed for patch allocation.
592 !> ... for the vertices
593 !! @ingroup fortran_interface
594 !!
595 !! NB: Avoid the case nblks=0 for empty patches, this might cause troubles
596 !! if a empty patch is used somewhere (and npromz gets wrong in the formulas below).
597 !!
598 !! from `mo_setup_subdivision`; name of ICON routine: `npromz_v`.
599 INTEGER(c_int) FUNCTION comin_descrdata_get_vert_npromz(jg) BIND(C)
600 INTEGER(c_int), INTENT(IN), VALUE :: jg ! domain index for comin_domain
601
602 comin_descrdata_get_vert_npromz = state%comin_descrdata_domain_data(jg)%verts%nverts - &
603 & (state%comin_descrdata_domain_data(jg)%verts%nblks-1)*state%comin_descrdata_global_data%nproma
605
606 !> Check if index is in bounds for an array pointer.
607 !! @internal
608 LOGICAL FUNCTION is_in_bounds_int(idx, arr) RESULT(is_in_bounds)
609 INTEGER, INTENT(IN) :: idx !< Array index.
610 INTEGER, POINTER, INTENT(IN) :: arr(:) !< Array. Must be a pointer to preserve bounds.
611
612 is_in_bounds = lbound(arr, dim=1) <= idx .AND. idx <= ubound(arr, dim=1)
613 END FUNCTION is_in_bounds_int
614
615 !> Check if index is in bounds for an array pointer.
616 !! @internal
617 LOGICAL FUNCTION is_in_bounds_ddd(idx, arr) RESULT(is_in_bounds)
618 INTEGER, INTENT(IN) :: idx !< Array index.
619 TYPE(t_comin_descrdata_domain_data), ALLOCATABLE, INTENT(IN) :: arr(:)
620 !^ Array. Must be a pointer/allocatable to preserve bounds.
621
622 is_in_bounds = lbound(arr, dim=1) <= idx .AND. idx <= ubound(arr, dim=1)
623 END FUNCTION is_in_bounds_ddd
624
625 SUBROUTINE comin_descrdata_to_c_signed_char_1(arr, ptr, arr_size)
626 INTEGER(c_signed_char), POINTER, INTENT(IN) :: arr(:)
627 TYPE(c_ptr), INTENT(OUT) :: ptr
628 INTEGER(C_INT), INTENT(OUT) :: arr_size(1)
629
630 ptr = c_loc(arr)
631 arr_size = shape(arr)
632 END SUBROUTINE comin_descrdata_to_c_signed_char_1
633
634 SUBROUTINE comin_descrdata_to_c_string(arr, ptr, arr_size)
635 CHARACTER(LEN=:), POINTER, INTENT(IN) :: arr
636 TYPE(c_ptr), INTENT(OUT) :: ptr
637 INTEGER(C_INT), INTENT(OUT) :: arr_size(1)
638
639 ptr = c_loc(arr)
640 arr_size = len_trim(arr)
641 END SUBROUTINE comin_descrdata_to_c_string
642
643 SUBROUTINE comin_descrdata_to_c_int_1(arr, ptr, arr_size)
644 INTEGER(C_INT), POINTER, INTENT(IN) :: arr(:)
645 TYPE(c_ptr), INTENT(OUT) :: ptr
646 INTEGER(C_INT), INTENT(OUT) :: arr_size(1)
647
648 ptr = c_loc(arr)
649 arr_size = shape(arr)
650 END SUBROUTINE comin_descrdata_to_c_int_1
651
652 SUBROUTINE comin_descrdata_to_c_int_2(arr, ptr, arr_size)
653 INTEGER(C_INT), POINTER, INTENT(IN) :: arr(:,:)
654 TYPE(c_ptr), INTENT(OUT) :: ptr
655 INTEGER(C_INT), INTENT(OUT) :: arr_size(2)
656
657 ptr = c_loc(arr)
658 arr_size = shape(arr)
659 END SUBROUTINE comin_descrdata_to_c_int_2
660
661 SUBROUTINE comin_descrdata_to_c_int_3(arr, ptr, arr_size)
662 INTEGER(C_INT), POINTER, INTENT(IN) :: arr(:,:,:)
663 TYPE(c_ptr), INTENT(OUT) :: ptr
664 INTEGER(C_INT), INTENT(OUT) :: arr_size(3)
665
666 ptr = c_loc(arr)
667 arr_size = shape(arr)
668 END SUBROUTINE comin_descrdata_to_c_int_3
669
670 SUBROUTINE comin_descrdata_to_c_real_1(arr, ptr, arr_size)
671 REAL(wp), POINTER, INTENT(IN) :: arr(:)
672 TYPE(c_ptr), INTENT(OUT) :: ptr
673 INTEGER(C_INT), INTENT(OUT) :: arr_size(1)
674
675 ptr = c_loc(arr)
676 arr_size = shape(arr)
677 END SUBROUTINE comin_descrdata_to_c_real_1
678
679 SUBROUTINE comin_descrdata_to_c_real_2(arr, ptr, arr_size)
680 REAL(wp), POINTER, INTENT(IN) :: arr(:,:)
681 TYPE(c_ptr), INTENT(OUT) :: ptr
682 INTEGER(C_INT), INTENT(OUT) :: arr_size(2)
683
684 ptr = c_loc(arr)
685 arr_size = shape(arr)
686 END SUBROUTINE comin_descrdata_to_c_real_2
687
688 SUBROUTINE comin_descrdata_to_c_real_3(arr, ptr, arr_size)
689 REAL(wp), POINTER, INTENT(IN) :: arr(:,:,:)
690 TYPE(c_ptr), INTENT(OUT) :: ptr
691 INTEGER(C_INT), INTENT(OUT) :: arr_size(3)
692
693 ptr = c_loc(arr)
694 arr_size = shape(arr)
695 END SUBROUTINE comin_descrdata_to_c_real_3
696
697 SUBROUTINE comin_descrdata_to_c_real_4(arr, ptr, arr_size)
698 REAL(wp), POINTER, INTENT(IN) :: arr(:,:,:,:)
699 TYPE(c_ptr), INTENT(OUT) :: ptr
700 INTEGER(C_INT), INTENT(OUT) :: arr_size(4)
701
702 ptr = c_loc(arr)
703 arr_size = shape(arr)
704 END SUBROUTINE comin_descrdata_to_c_real_4
705
706END MODULE comin_descrdata
void comin_error_set(t_comin_error_code error_code)
void comin_state_set_simulation_interval(const char *exp_start, const char *exp_stop, const char *run_start, const char *run_stop)
const char * comin_descrdata_get_simulation_interval_exp_start()
const char * comin_descrdata_get_simulation_interval_exp_stop()
const char * comin_descrdata_get_simulation_interval_run_start()
const char * comin_descrdata_get_simulation_interval_run_stop()
Global data is invariant wrt the computational grid and never changed or updated.
Conversion of global edge index to MPI-process local index.
Patch grid data structure, gathering information on grids.
Simulation status information, sim_current contains current time step.
Conversion of global vert index to MPI-process local index.
Conversion of global cell index to MPI-process local index.
Receive pointer on array storing timestep information for all domains.
subroutine, public comin_descrdata_get_vert_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk)
Computes the start and end block indices for loops over vertex-based variables.
integer(c_int) function, public comin_descrdata_get_vert_npromz(jg)
Calculate npromz value for the blocking, needed for patch allocation. ... for the vertices.
subroutine, public comin_descrdata_get_vert_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, i_endidx, irl_start, irl_end)
Computes the start and end indices of do loops for vertex-based variables.
subroutine, public comin_descrdata_get_cell_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, i_endidx, irl_start, irl_end)
Computes the start and end indices of do loops for cell-based variables.
subroutine, public comin_descrdata_get_edge_indices(jg, i_blk, i_startblk, i_endblk, i_startidx, i_endidx, irl_start, irl_end)
Computes the start and end indices of do loops for edge-based variables.
integer(c_int) function, public comin_descrdata_get_cell_npromz(jg)
Calculate npromz value for the blocking, needed for patch allocation. ... for the cells.
integer, parameter wp
working precision
integer(c_int) function, public comin_descrdata_get_index(idx1d)
Auxiliary function: conversion of 1D to 2D indices.
subroutine, public comin_descrdata_get_edge_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk)
Computes the start and end block indices for loops over edge-based variables.
integer(c_int) function, public comin_descrdata_get_edge_npromz(jg)
Calculate npromz value for the blocking, needed for patch allocation. ... for the edges.
subroutine, public comin_descrdata_get_cell_block_limits(jg, irl_start, irl_end, i_startblk, i_endblk)
Computes the start and end block indices for loops over cell-based variables.
type(t_comin_descrdata_simulation_interval) function, public comin_descrdata_get_simulation_interval()
request a pointer to simulation status
integer(c_int) function, public comin_descrdata_get_block(idx1d)
auxiliary functions taken from ICON, version 2.6.5
subroutine, public comin_current_get_datetime(sim_time_current)
Retrieve time stamp info, current time information.
subroutine, public comin_descrdata_set_simulation_interval(exp_start, exp_stop, run_start, run_stop)
Fill time stamp info.
subroutine, public comin_descrdata_finalize()
Clean descriptive data structure in ComIn currently no content but keep for future use.
subroutine, public comin_descrdata_set_domain(comin_domain_info)
Set up data type for grid data.
subroutine, public comin_current_set_datetime(sim_time_current)
Update time stamp info, current time information.
subroutine, public comin_descrdata_set_global(comin_global_info)
Fill global data.
type(t_comin_state), pointer, public state