YAC 3.7.1
Yet Another Coupler
Loading...
Searching...
No Matches
yac_module.F90
Go to the documentation of this file.
1! Copyright (c) 2024 The YAC Authors
2!
3! SPDX-License-Identifier: BSD-3-Clause
4
5#ifdef HAVE_CONFIG_H
6! Get the definition of the 'YAC_MPI_FINT_FC_KIND' macro.
7#include "config.h"
8#endif
9
10module yac
11 use, intrinsic :: iso_c_binding, only : c_int, c_long, &
12 & c_long_long, c_short, c_char
13
14 !----------------------------------------------------------------------
18 !----------------------------------------------------------------------
19
20 public
21
22 integer, parameter :: yac_max_charlen = 132
23 integer, parameter :: yac_mpi_fint_kind = yac_mpi_fint_fc_kind
24 double precision, parameter :: yac_frac_mask_no_value = 133713371337.0d0
25
26 !----------------------------------------------------------------------
30 !----------------------------------------------------------------------
31
32 enum, bind(c)
33 enumerator :: yac_location_cell = 0
34 enumerator :: yac_location_corner = 1
35 enumerator :: yac_location_edge = 2
36 end enum
37
38 !----------------------------------------------------------------------
42 !----------------------------------------------------------------------
43
44 enum, bind(c)
45 enumerator :: yac_exchange_type_none = 0
46 enumerator :: yac_exchange_type_source = 1
47 enumerator :: yac_exchange_type_target = 2
48 end enum
49
50 !----------------------------------------------------------------------
54 !----------------------------------------------------------------------
55
56 enum, bind(c)
57 enumerator :: yac_action_none = 0
58 enumerator :: yac_action_reduction = 1
59 enumerator :: yac_action_coupling = 2
62 enumerator :: yac_action_out_of_bound = 8
63 end enum
64
65 !----------------------------------------------------------------------
69 !----------------------------------------------------------------------
70 enum, bind(c)
71 enumerator :: yac_time_unit_millisecond = 0
72 enumerator :: yac_time_unit_second = 1
73 enumerator :: yac_time_unit_minute = 2
74 enumerator :: yac_time_unit_hour = 3
75 enumerator :: yac_time_unit_day = 4
76 enumerator :: yac_time_unit_month = 5
77 enumerator :: yac_time_unit_year = 6
78 enumerator :: yac_time_unit_iso_format = 7
79 end enum
80
81 !----------------------------------------------------------------------
85 !----------------------------------------------------------------------
86
87 enum, bind(c)
88 enumerator :: yac_reduction_time_none = 0
93 end enum
94
95 !----------------------------------------------------------------------
99 !----------------------------------------------------------------------
100 enum, bind(c)
101 enumerator :: yac_calendar_not_set = 0
102 enumerator :: yac_proleptic_gregorian = 1
103 enumerator :: yac_year_of_365_days = 2
104 enumerator :: yac_year_of_360_days = 3
105 end enum
106
107 !----------------------------------------------------------------------
111 !----------------------------------------------------------------------
112
113 enum, bind(c)
114 enumerator :: yac_avg_arithmetic = 0
115 enumerator :: yac_avg_dist = 1
116 enumerator :: yac_avg_bary = 2
117 end enum
118
119 enum, bind(c)
120 enumerator :: yac_ncc_avg = 0
121 enumerator :: yac_ncc_dist = 1
122 end enum
123
124 enum, bind(c)
125 enumerator :: yac_nnn_avg = 0
126 enumerator :: yac_nnn_dist = 1
127 enumerator :: yac_nnn_gauss = 2
128 enumerator :: yac_nnn_rbf = 3
129 end enum
130
131 enum, bind(c)
132 enumerator :: yac_conserv_destarea = 0
133 enumerator :: yac_conserv_fracarea = 1
134 end enum
135
136 enum, bind(c)
137 enumerator :: yac_spmap_avg = 0
138 enumerator :: yac_spmap_dist = 1
139 end enum
140
141 enum, bind(c)
142 enumerator :: yac_spmap_none = 0
143 enumerator :: yac_spmap_srcarea = 1
144 enumerator :: yac_spmap_invtgtarea = 2
145 enumerator :: yac_spmap_fracarea = 3
146 end enum
147
148 enum, bind(c)
149 enumerator :: yac_file_missing_error = 0
150 enumerator :: yac_file_missing_cont = 1
151 end enum
152
153 enum, bind(c)
154 enumerator :: yac_file_success_stop = 0
157 enumerator :: yac_file_success_cont = 1
159 end enum
160
161 !----------------------------------------------------------------------
165 !----------------------------------------------------------------------
166
167 enum, bind(c)
170 end enum
171
172 enum, bind(c)
176 end enum
177
178 !----------------------------------
179 ! Handling of existing weight files
180 !----------------------------------
181
182 enum, bind(c)
184 enumerator :: yac_wgt_on_existing_keep = 1
186 end enum
187
188 !----------------------------------------------------------------------
192 !----------------------------------------------------------------------
193
196
197 !----------------------------------------------------------------------
201 !----------------------------------------------------------------------
202
204 real, pointer :: p(:)
205 end type yac_real_ptr
206
208 double precision, pointer :: p(:)
209 end type yac_dble_ptr
210
212 character(len=:), allocatable :: string
213 end type yac_string
214
215 !----------------------------------------------------------------------
219 !----------------------------------------------------------------------
220
222 subroutine yac_fmpi_handshake ( comm, group_names, group_comms )
223 import :: yac_max_charlen
224 integer, intent(in) :: comm
225 character(len=YAC_MAX_CHARLEN), intent(in) :: &
226 group_names(:)
227 integer, intent(out) :: &
228 group_comms(SIZE(group_names))
229 end subroutine yac_fmpi_handshake
230 end interface yac_fmpi_handshake
231
232 !----------------------------------------------------------------------
236 !----------------------------------------------------------------------
237
239
240 subroutine yac_finit_comm ( comm )
241
242 integer, intent(in) :: comm
243
244 end subroutine yac_finit_comm
245
246 subroutine yac_finit_comm_instance ( comm, yac_instance_id )
247
248 integer, intent(in) :: comm
249 integer, intent(out) :: yac_instance_id
250
251 end subroutine yac_finit_comm_instance
252
253 end interface yac_finit_comm
254
255 interface yac_finit
256
257 subroutine yac_finit ( )
258
259 end subroutine yac_finit
260
261 subroutine yac_finit_instance ( yac_instance_id)
262
263 integer, intent(out) :: yac_instance_id
264
265 end subroutine yac_finit_instance
266
267 end interface yac_finit
268
270 subroutine yac_finit_dummy ()
271 end subroutine yac_finit_dummy
272
273 subroutine yac_finit_comm_dummy ( comm )
274
275 integer, intent(in) :: comm
276
277 end subroutine yac_finit_comm_dummy
278 end interface yac_finit_dummy
279
280 !----------------------------------------------------------------------
284 !----------------------------------------------------------------------
285
286 interface
292 end interface
293
294 !----------------------------------------------------------------------
298 !----------------------------------------------------------------------
299
301 subroutine yac_fread_config_yaml (yaml_filename)
302 character(len=*), intent(in) :: yaml_filename
303 end subroutine yac_fread_config_yaml
304
305 subroutine yac_fread_config_yaml_instance (yac_instance_id, yaml_filename)
306 integer, intent(in) :: yac_instance_id
307 character(len=*), intent(in) :: yaml_filename
308 end subroutine yac_fread_config_yaml_instance
309 end interface yac_fread_config_yaml
310
312 subroutine yac_fread_config_json (json_filename)
313 character(len=*), intent(in) :: json_filename
314 end subroutine yac_fread_config_json
315
316 subroutine yac_fread_config_json_instance (yac_instance_id, json_filename)
317 integer, intent(in) :: yac_instance_id
318 character(len=*), intent(in) :: json_filename
319 end subroutine yac_fread_config_json_instance
320 end interface yac_fread_config_json
321
322 !----------------------------------------------------------------------
326 !----------------------------------------------------------------------
327
330 filename, fileformat, sync_location, include_definitions)
331 character(len=*), intent(in) :: filename
332 integer, intent(in) :: fileformat
333 integer, intent(in) :: sync_location
334 logical, intent(in), optional :: include_definitions
336 end subroutine yac_fset_config_output_file
337
339 yac_instance_id, filename, fileformat, sync_location, &
340 include_definitions)
341 integer, intent(in) :: yac_instance_id
342 character(len=*), intent(in) :: filename
343 integer, intent(in) :: fileformat
344 integer, intent(in) :: sync_location
345 logical, intent(in), optional :: include_definitions
348 end interface yac_fset_config_output_file
349
350 !----------------------------------------------------------------------
354 !----------------------------------------------------------------------
355
358 gridname, filename)
359 character(len=*), intent(in) :: gridname
360 character(len=*), intent(in) :: filename
361 end subroutine yac_fset_grid_output_file
362
364 yac_instance_id, gridname, filename)
365 integer, intent(in) :: yac_instance_id
366 character(len=*), intent(in) :: gridname
367 character(len=*), intent(in) :: filename
369 end interface yac_fset_grid_output_file
370
371 !----------------------------------------------------------------------
375 !----------------------------------------------------------------------
376
377 interface yac_fcleanup
378
379 subroutine yac_fcleanup ()
380 end subroutine yac_fcleanup
381
382 subroutine yac_fcleanup_instance ( yac_instance_id )
383
384 integer, intent(in) :: yac_instance_id
385
386 end subroutine yac_fcleanup_instance
387
388 end interface yac_fcleanup
389
390 !----------------------------------------------------------------------
394 !----------------------------------------------------------------------
395
397
398 subroutine yac_ffinalize ()
399 end subroutine yac_ffinalize
400
401 subroutine yac_ffinalize_instance ( yac_instance_id )
402
403 integer, intent(in) :: yac_instance_id
404
405 end subroutine
406
407 end interface yac_ffinalize
408
409 !----------------------------------------------------------------------
413 !----------------------------------------------------------------------
414
415 interface
416 function yac_fget_version () result (version_string)
417 character(len=:), ALLOCATABLE :: version_string
418 end function yac_fget_version
419 end interface
420
421 !----------------------------------------------------------------------
425 !----------------------------------------------------------------------
426
428
429 SUBROUTINE yac_fpredef_comp ( comp_name, comp_id )
430
431 character(len=*), intent(in) :: comp_name
432 integer, intent(out) :: comp_id
433
434 END SUBROUTINE yac_fpredef_comp
435
436 SUBROUTINE yac_fpredef_comp_instance ( yac_instance_id, comp_name, comp_id )
437
438 integer, intent(in) :: yac_instance_id
439 character(len=*), intent(in) :: comp_name
440 integer, intent(out) :: comp_id
441
442 END SUBROUTINE yac_fpredef_comp_instance
443
444 END INTERFACE yac_fpredef_comp
445
447
448 subroutine yac_fdef_comp ( comp_name, comp_id )
449
450 character(len=*), intent(in) :: comp_name
451 integer, intent(out) :: comp_id
452
453 end subroutine yac_fdef_comp
454
455 subroutine yac_fdef_comp_instance ( yac_instance_id, comp_name, comp_id )
456
457 integer, intent(in) :: yac_instance_id
458 character(len=*), intent(in) :: comp_name
459 integer, intent(out) :: comp_id
460
461 end subroutine yac_fdef_comp_instance
462
463 end interface yac_fdef_comp
464
466
467 subroutine yac_fdef_comps ( comp_names, num_comps, comp_ids )
468
469 use, intrinsic :: iso_c_binding, only : c_char
470 import :: yac_max_charlen
471
472 integer, intent(in) :: num_comps
473 character(kind=c_char, len=YAC_MAX_CHARLEN), intent(in) :: &
474 comp_names(num_comps)
475 integer, intent(out) :: comp_ids(num_comps)
476
477 end subroutine yac_fdef_comps
478
479 subroutine yac_fdef_comps_instance ( yac_instance_id, &
480 comp_names, &
481 num_comps, &
482 comp_ids )
483
484 use, intrinsic :: iso_c_binding, only : c_char
485 import :: yac_max_charlen
486
487 integer, intent(in) :: yac_instance_id
488 integer, intent(in) :: num_comps
489 character(kind=c_char, len=YAC_MAX_CHARLEN), intent(in) :: &
490 comp_names(num_comps)
491 integer, intent(out) :: comp_ids(num_comps)
492
493 end subroutine yac_fdef_comps_instance
494
495 end interface yac_fdef_comps
496
498
499 subroutine yac_fdef_comp_dummy ( )
500 end subroutine yac_fdef_comp_dummy
501
502 subroutine yac_fdef_comp_dummy_instance ( yac_instance_id )
503
504 integer, intent(in) :: yac_instance_id
505
506 end subroutine yac_fdef_comp_dummy_instance
507
508 end interface yac_fdef_comp_dummy
509
510 !----------------------------------------------------------------------
514 !----------------------------------------------------------------------
515
517
518 subroutine yac_fdef_datetime ( start_datetime, end_datetime )
519
520 character(len=*), intent(in), optional :: start_datetime
521 character(len=*), intent(in), optional :: end_datetime
522
523 end subroutine yac_fdef_datetime
524
525 subroutine yac_fdef_datetime_instance ( yac_instance_id, &
526 start_datetime, &
527 end_datetime )
528
529 integer, intent(in) :: yac_instance_id
530 character(len=*), intent(in), optional :: start_datetime
531 character(len=*), intent(in), optional :: end_datetime
532
533 end subroutine yac_fdef_datetime_instance
534
535 end interface yac_fdef_datetime
536
538
539 subroutine yac_fdef_calendar( calendar )
540 integer, intent(in) :: calendar
541 end subroutine yac_fdef_calendar
542
543 end interface yac_fdef_calendar
544
546
547 subroutine yac_fget_calendar( calendar )
548 integer, intent(out) :: calendar
549 end subroutine yac_fget_calendar
550
551 end interface yac_fget_calendar
552
553 !----------------------------------------------------------------------
557 !----------------------------------------------------------------------
558
560
561 function yac_fget_start_datetime () result (start_datetime_string)
562 character(len=:), ALLOCATABLE :: start_datetime_string
563 end function yac_fget_start_datetime
564
565 function yac_fget_start_datetime_instance (yac_instance_id) &
566 result(start_datetime_string)
567 integer, intent(in) :: yac_instance_id
568 character(len=:), ALLOCATABLE :: start_datetime_string
570
571 end interface yac_fget_start_datetime
572
574
575 function yac_fget_end_datetime () result (end_datetime_string)
576 character(len=:), ALLOCATABLE :: end_datetime_string
577 end function yac_fget_end_datetime
578
579 function yac_fget_end_datetime_instance (yac_instance_id) &
580 result(end_datetime_string)
581 integer, intent(in) :: yac_instance_id
582 character(len=:), ALLOCATABLE :: end_datetime_string
584
585 end interface yac_fget_end_datetime
586
587 !----------------------------------------------------------------------
591 !----------------------------------------------------------------------
592
594
614 subroutine yac_fdef_grid_nonuniform_real ( grid_name, &
615 nbr_vertices, &
616 nbr_cells, &
617 nbr_connections, &
618 nbr_vertices_per_cell, &
619 x_vertices, &
620 y_vertices, &
621 cell_to_vertex, &
622 grid_id, &
623 use_ll_edges)
624
625 character(len=*), intent(in) :: grid_name
626 integer, intent(in) :: nbr_vertices
627 integer, intent(in) :: nbr_cells
628 integer, intent(in) :: nbr_connections
629 integer, intent(in) :: nbr_vertices_per_cell(nbr_cells)
630 real, intent(in) :: x_vertices(nbr_vertices)
631 real, intent(in) :: y_vertices(nbr_vertices)
632 integer, intent(in) :: cell_to_vertex(nbr_connections)
633 integer, intent(out) :: grid_id
634 logical, optional, intent(in) :: use_ll_edges
635
636 end subroutine yac_fdef_grid_nonuniform_real
637
657 subroutine yac_fdef_grid_nonuniform_dble ( grid_name, &
658 nbr_vertices, &
659 nbr_cells, &
660 nbr_connections, &
661 nbr_vertices_per_cell, &
662 x_vertices, &
663 y_vertices, &
664 cell_to_vertex, &
665 grid_id, &
666 use_ll_edges)
667
668 character(len=*), intent(in) :: grid_name
669 integer, intent(in) :: nbr_vertices
670 integer, intent(in) :: nbr_cells
671 integer, intent(in) :: nbr_connections
672 integer, intent(in) :: nbr_vertices_per_cell(nbr_cells)
673 double precision, intent(in) :: x_vertices(nbr_vertices)
674 double precision, intent(in) :: y_vertices(nbr_vertices)
675 integer, intent(in) :: cell_to_vertex(nbr_connections)
676 integer, intent(out) :: grid_id
677 logical, optional, intent(in) :: use_ll_edges
678
679 end subroutine yac_fdef_grid_nonuniform_dble
680
699 subroutine yac_fdef_grid_unstruct_real ( grid_name, &
700 nbr_vertices, &
701 nbr_cells, &
702 nbr_vertices_per_cell, &
703 x_vertices, &
704 y_vertices, &
705 cell_to_vertex, &
706 grid_id, &
707 use_ll_edges)
708
709 character(len=*), intent(in) :: grid_name
710 integer, intent(in) :: nbr_vertices
711 integer, intent(in) :: nbr_cells
712 integer, intent(in) :: nbr_vertices_per_cell
713 real, intent(in) :: x_vertices(nbr_vertices)
714 real, intent(in) :: y_vertices(nbr_vertices)
715 integer, intent(in) :: cell_to_vertex( &
716 nbr_vertices_per_cell,nbr_cells)
717 integer, intent(out) :: grid_id
718 logical, optional, intent(in) :: use_ll_edges
719
720 end subroutine yac_fdef_grid_unstruct_real
721
740 subroutine yac_fdef_grid_unstruct_dble ( grid_name, &
741 nbr_vertices, &
742 nbr_cells, &
743 nbr_vertices_per_cell, &
744 x_vertices, &
745 y_vertices, &
746 cell_to_vertex, &
747 grid_id, &
748 use_ll_edges)
749
750 character(len=*), intent(in) :: grid_name
751 integer, intent(in) :: nbr_vertices
752 integer, intent(in) :: nbr_cells
753 integer, intent(in) :: nbr_vertices_per_cell
754 double precision, intent(in) :: x_vertices(nbr_vertices)
755 double precision, intent(in) :: y_vertices(nbr_vertices)
756 integer, intent(in) :: cell_to_vertex( &
757 nbr_vertices_per_cell,nbr_cells)
758 integer, intent(out) :: grid_id
759 logical, optional, intent(in) :: use_ll_edges
760
761 end subroutine yac_fdef_grid_unstruct_dble
762
784 subroutine yac_fdef_grid_nonuniform_edge_real ( grid_name, &
785 nbr_vertices, &
786 nbr_cells, &
787 nbr_edges, &
788 nbr_connections, &
789 nbr_edges_per_cell, &
790 x_vertices, &
791 y_vertices, &
792 cell_to_edge, &
793 edge_to_vertex, &
794 grid_id, &
795 use_ll_edges)
796
797 character(len=*), intent(in) :: grid_name
798 integer, intent(in) :: nbr_vertices
799 integer, intent(in) :: nbr_cells
800 integer, intent(in) :: nbr_edges
801 integer, intent(in) :: nbr_connections
802 integer, intent(in) :: nbr_edges_per_cell(nbr_cells)
803 real, intent(in) :: x_vertices(nbr_vertices)
804 real, intent(in) :: y_vertices(nbr_vertices)
805 integer, intent(in) :: cell_to_edge(nbr_connections)
806 integer, intent(in) :: edge_to_vertex(2,nbr_edges)
807 integer, intent(out) :: grid_id
808 logical, optional, intent(in) :: use_ll_edges
809
811
833 subroutine yac_fdef_grid_nonuniform_edge_dble ( grid_name, &
834 nbr_vertices, &
835 nbr_cells, &
836 nbr_edges, &
837 nbr_connections, &
838 nbr_edges_per_cell, &
839 x_vertices, &
840 y_vertices, &
841 cell_to_edge, &
842 edge_to_vertex, &
843 grid_id, &
844 use_ll_edges)
845
846 character(len=*), intent(in) :: grid_name
847 integer, intent(in) :: nbr_vertices
848 integer, intent(in) :: nbr_cells
849 integer, intent(in) :: nbr_edges
850 integer, intent(in) :: nbr_connections
851 integer, intent(in) :: nbr_edges_per_cell(nbr_cells)
852 double precision, intent(in) :: x_vertices(nbr_vertices)
853 double precision, intent(in) :: y_vertices(nbr_vertices)
854 integer, intent(in) :: cell_to_edge(nbr_connections)
855 integer, intent(in) :: edge_to_vertex(2,nbr_edges)
856 integer, intent(out) :: grid_id
857 logical, optional, intent(in) :: use_ll_edges
858
860
881 subroutine yac_fdef_grid_unstruct_edge_real ( grid_name, &
882 nbr_vertices, &
883 nbr_cells, &
884 nbr_edges, &
885 nbr_edges_per_cell, &
886 x_vertices, &
887 y_vertices, &
888 cell_to_edge, &
889 edge_to_vertex, &
890 grid_id, &
891 use_ll_edges)
892
893 character(len=*), intent(in) :: grid_name
894 integer, intent(in) :: nbr_vertices
895 integer, intent(in) :: nbr_cells
896 integer, intent(in) :: nbr_edges
897 integer, intent(in) :: nbr_edges_per_cell
898 real, intent(in) :: x_vertices(nbr_vertices)
899 real, intent(in) :: y_vertices(nbr_vertices)
900 integer, intent(in) :: cell_to_edge( &
901 nbr_edges_per_cell,nbr_cells)
902 integer, intent(in) :: edge_to_vertex(2,nbr_edges)
903 integer, intent(out) :: grid_id
904 logical, optional, intent(in) :: use_ll_edges
905
907
928 subroutine yac_fdef_grid_unstruct_edge_dble ( grid_name, &
929 nbr_vertices, &
930 nbr_cells, &
931 nbr_edges, &
932 nbr_edges_per_cell, &
933 x_vertices, &
934 y_vertices, &
935 cell_to_edge, &
936 edge_to_vertex, &
937 grid_id, &
938 use_ll_edges)
939
940 character(len=*), intent(in) :: grid_name
941 integer, intent(in) :: nbr_vertices
942 integer, intent(in) :: nbr_cells
943 integer, intent(in) :: nbr_edges
944 integer, intent(in) :: nbr_edges_per_cell
945 double precision, intent(in) :: x_vertices(nbr_vertices)
946 double precision, intent(in) :: y_vertices(nbr_vertices)
947 integer, intent(in) :: cell_to_edge( &
948 nbr_edges_per_cell,nbr_cells)
949 integer, intent(in) :: edge_to_vertex(2,nbr_edges)
950 integer, intent(out) :: grid_id
951 logical, optional, intent(in) :: use_ll_edges
952
954
962 subroutine yac_fdef_grid_curve2d_real ( grid_name, &
963 nbr_vertices, &
964 cyclic, &
965 x_vertices, &
966 y_vertices, &
967 grid_id)
968
969 character(len=*), intent(in) :: grid_name
970 integer, intent(in) :: nbr_vertices(2)
971 integer, intent(in) :: cyclic(2)
972 real, intent(in) :: &
973 x_vertices(nbr_vertices(1),nbr_vertices(2))
974 real, intent(in) :: &
975 y_vertices(nbr_vertices(1),nbr_vertices(2))
976 integer, intent(out) :: grid_id
977
978 end subroutine yac_fdef_grid_curve2d_real
979
987 subroutine yac_fdef_grid_curve2d_dble ( grid_name, &
988 nbr_vertices, &
989 cyclic, &
990 x_vertices, &
991 y_vertices, &
992 grid_id)
993
994 character(len=*), intent(in) :: grid_name
995 integer, intent(in) :: nbr_vertices(2)
996 integer, intent(in) :: cyclic(2)
997 double precision, intent(in) :: &
998 x_vertices(nbr_vertices(1),nbr_vertices(2))
999 double precision, intent(in) :: &
1000 y_vertices(nbr_vertices(1),nbr_vertices(2))
1001 integer, intent(out) :: grid_id
1002
1003 end subroutine yac_fdef_grid_curve2d_dble
1004
1012 subroutine yac_fdef_grid_reg2d_real ( grid_name, &
1013 nbr_vertices, &
1014 cyclic, &
1015 x_vertices, &
1016 y_vertices, &
1017 grid_id)
1018
1019 character(len=*), intent(in) :: grid_name
1020 integer, intent(in) :: nbr_vertices(2)
1021 integer, intent(in) :: cyclic(2)
1022 real, intent(in) :: x_vertices(nbr_vertices(1))
1023 real, intent(in) :: y_vertices(nbr_vertices(2))
1024 integer, intent(out) :: grid_id
1025
1026 end subroutine yac_fdef_grid_reg2d_real
1027
1035 subroutine yac_fdef_grid_reg2d_dble ( grid_name, &
1036 nbr_vertices, &
1037 cyclic, &
1038 x_vertices, &
1039 y_vertices, &
1040 grid_id)
1041
1042 character(len=*), intent(in) :: grid_name
1043 integer, intent(in) :: nbr_vertices(2)
1044 integer, intent(in) :: cyclic(2)
1045 double precision, intent(in) :: x_vertices(nbr_vertices(1))
1046 double precision, intent(in) :: y_vertices(nbr_vertices(2))
1047 integer, intent(out) :: grid_id
1048
1049 end subroutine yac_fdef_grid_reg2d_dble
1050
1057 subroutine yac_fdef_grid_cloud_real ( grid_name, &
1058 nbr_points, &
1059 x_points, &
1060 y_points, &
1061 grid_id)
1062
1063 character(len=*), intent(in) :: grid_name
1064 integer, intent(in) :: nbr_points
1065 real, intent(in) :: x_points(nbr_points)
1066 real, intent(in) :: y_points(nbr_points)
1067 integer, intent(out) :: grid_id
1068
1069 end subroutine yac_fdef_grid_cloud_real
1070
1077 subroutine yac_fdef_grid_cloud_dble ( grid_name, &
1078 nbr_points, &
1079 x_points, &
1080 y_points, &
1081 grid_id)
1082
1083 character(len=*), intent(in) :: grid_name
1084 integer, intent(in) :: nbr_points
1085 double precision, intent(in) :: x_points(nbr_points)
1086 double precision, intent(in) :: y_points(nbr_points)
1087 integer, intent(out) :: grid_id
1088
1089 end subroutine yac_fdef_grid_cloud_dble
1090
1100 subroutine yac_fdef_grid_reg2d_rot_real ( grid_name, &
1101 nbr_vertices, &
1102 cyclic, &
1103 x_vertices, &
1104 y_vertices, &
1105 x_north_pole, &
1106 y_north_pole, &
1107 grid_id)
1108
1109 character(len=*), intent(in) :: grid_name
1110 integer, intent(in) :: nbr_vertices(2)
1111 integer, intent(in) :: cyclic(2)
1112 real, intent(in) :: x_vertices(nbr_vertices(1))
1113 real, intent(in) :: y_vertices(nbr_vertices(2))
1114 real, intent(in) :: x_north_pole
1115 real, intent(in) :: y_north_pole
1116 integer, intent(out) :: grid_id
1117
1118 end subroutine yac_fdef_grid_reg2d_rot_real
1119
1129 subroutine yac_fdef_grid_reg2d_rot_dble ( grid_name, &
1130 nbr_vertices, &
1131 cyclic, &
1132 x_vertices, &
1133 y_vertices, &
1134 x_north_pole, &
1135 y_north_pole, &
1136 grid_id)
1137
1138 character(len=*), intent(in) :: grid_name
1139 integer, intent(in) :: nbr_vertices(2)
1140 integer, intent(in) :: cyclic(2)
1141 double precision, intent(in) :: x_vertices(nbr_vertices(1))
1142 double precision, intent(in) :: y_vertices(nbr_vertices(2))
1143 double precision, intent(in) :: x_north_pole
1144 double precision, intent(in) :: y_north_pole
1145 integer, intent(out) :: grid_id
1146
1147 end subroutine yac_fdef_grid_reg2d_rot_dble
1148
1149 end interface yac_fdef_grid
1150
1151 !----------------------------------------------------------------------
1155 !----------------------------------------------------------------------
1156
1158
1159 subroutine yac_fdef_points_reg2d_real ( grid_id, &
1160 nbr_points, &
1161 location, &
1162 x_points_real, &
1163 y_points_real, &
1164 point_id )
1165
1166 integer, intent(in) :: grid_id
1167 integer, intent(in) :: nbr_points(2)
1168 integer, intent(in) :: location
1169 real, intent(in) :: x_points_real(nbr_points(1))
1170 real, intent(in) :: y_points_real(nbr_points(2))
1171 integer, intent(out) :: point_id
1172
1173 end subroutine yac_fdef_points_reg2d_real
1174
1175 subroutine yac_fdef_points_reg2d_dble ( grid_id, &
1176 nbr_points, &
1177 location, &
1178 x_points, &
1179 y_points, &
1180 point_id )
1181
1182 integer, intent(in) :: grid_id
1183 integer, intent(in) :: nbr_points(2)
1184 integer, intent(in) :: location
1185 double precision, intent(in) :: x_points(nbr_points(1))
1186 double precision, intent(in) :: y_points(nbr_points(2))
1187 integer, intent(out) :: point_id
1188
1189 end subroutine yac_fdef_points_reg2d_dble
1190
1191 subroutine yac_fdef_points_curve2d_real ( grid_id, &
1192 nbr_points, &
1193 location, &
1194 x_points_real, &
1195 y_points_real, &
1196 point_id )
1197
1198 integer, intent(in) :: grid_id
1199 integer, intent(in) :: nbr_points(2)
1200 integer, intent(in) :: location
1201 real, intent(in) :: &
1202 x_points_real(nbr_points(1),nbr_points(2))
1203 real, intent(in) :: &
1204 y_points_real(nbr_points(1),nbr_points(2))
1205 integer, intent(out) :: point_id
1206
1207 end subroutine yac_fdef_points_curve2d_real
1208
1209 subroutine yac_fdef_points_curve2d_dble ( grid_id, &
1210 nbr_points, &
1211 location, &
1212 x_points, &
1213 y_points, &
1214 point_id )
1215
1216 integer, intent(in) :: grid_id
1217 integer, intent(in) :: nbr_points(2)
1218 integer, intent(in) :: location
1219 double precision, intent(in) :: &
1220 x_points(nbr_points(1),nbr_points(2))
1221 double precision, intent(in) :: &
1222 y_points(nbr_points(1),nbr_points(2))
1223 integer, intent(out) :: point_id
1224
1225 end subroutine yac_fdef_points_curve2d_dble
1226
1227 subroutine yac_fdef_points_unstruct_real ( grid_id, &
1228 nbr_points, &
1229 location, &
1230 x_points_real, &
1231 y_points_real, &
1232 point_id )
1233
1234 integer, intent(in) :: grid_id
1235 integer, intent(in) :: nbr_points
1236 integer, intent(in) :: location
1237 real, intent(in) :: x_points_real(nbr_points)
1238 real, intent(in) :: y_points_real(nbr_points)
1239 integer, intent(out) :: point_id
1240
1241 end subroutine yac_fdef_points_unstruct_real
1242
1243 subroutine yac_fdef_points_unstruct_dble ( grid_id, &
1244 nbr_points, &
1245 location, &
1246 x_points, &
1247 y_points, &
1248 point_id )
1249
1250 integer, intent(in) :: grid_id
1251 integer, intent(in) :: nbr_points
1252 integer, intent(in) :: location
1253 double precision, intent(in) :: x_points(nbr_points)
1254 double precision, intent(in) :: y_points(nbr_points)
1255 integer, intent(out) :: point_id
1256
1257 end subroutine yac_fdef_points_unstruct_dble
1258
1259 subroutine yac_fdef_points_reg2d_rot_real ( grid_id, &
1260 nbr_points, &
1261 location, &
1262 x_points_real, &
1263 y_points_real, &
1264 x_north_pole_real, &
1265 y_north_pole_real, &
1266 point_id )
1267
1268 integer, intent(in) :: grid_id
1269 integer, intent(in) :: nbr_points(2)
1270 integer, intent(in) :: location
1271 real, intent(in) :: x_points_real(nbr_points(1))
1272 real, intent(in) :: y_points_real(nbr_points(2))
1273 real, intent(in) :: x_north_pole_real
1274 real, intent(in) :: y_north_pole_real
1275 integer, intent(out) :: point_id
1276
1277 end subroutine yac_fdef_points_reg2d_rot_real
1278
1279 subroutine yac_fdef_points_reg2d_rot_dble ( grid_id, &
1280 nbr_points, &
1281 location, &
1282 x_points, &
1283 y_points, &
1284 x_north_pole, &
1285 y_north_pole, &
1286 point_id )
1287
1288 integer, intent(in) :: grid_id
1289 integer, intent(in) :: nbr_points(2)
1290 integer, intent(in) :: location
1291 double precision, intent(in) :: x_points(nbr_points(1))
1292 double precision, intent(in) :: y_points(nbr_points(2))
1293 double precision, intent(in) :: x_north_pole
1294 double precision, intent(in) :: y_north_pole
1295 integer, intent(out) :: point_id
1296
1297 end subroutine yac_fdef_points_reg2d_rot_dble
1298
1299 end interface yac_fdef_points
1300
1301 !----------------------------------------------------------------------
1305 !----------------------------------------------------------------------
1306
1307 interface
1308
1309 subroutine yac_fset_global_index ( global_index, &
1310 location, &
1311 grid_id )
1312
1313 integer, intent(in) :: global_index(*)
1314 integer, intent(in) :: location
1315 integer, intent(in) :: grid_id
1316
1317 end subroutine yac_fset_global_index
1318
1319 end interface
1320
1321 !----------------------------------------------------------------------
1325 !----------------------------------------------------------------------
1326
1328
1329 subroutine yac_fset_core_imask ( is_core, &
1330 location, &
1331 grid_id )
1332
1333 integer, intent(in) :: is_core(*)
1334 integer, intent(in) :: location
1335 integer, intent(in) :: grid_id
1336
1337 end subroutine yac_fset_core_imask
1338
1339 subroutine yac_fset_core_lmask ( is_core, &
1340 location, &
1341 grid_id )
1342
1343 logical, intent(in) :: is_core(*)
1344 integer, intent(in) :: location
1345 integer, intent(in) :: grid_id
1346
1347 end subroutine yac_fset_core_lmask
1348
1349 end interface yac_fset_core_mask
1350
1351 !----------------------------------------------------------------------
1355 !----------------------------------------------------------------------
1356
1358
1359 subroutine yac_fset_imask ( is_valid, points_id )
1360
1361 integer, intent(in) :: is_valid(*)
1362 integer, intent(in) :: points_id
1363
1364 end subroutine yac_fset_imask
1365
1366 subroutine yac_fset_lmask ( is_valid, points_id )
1367
1368 logical, intent(in) :: is_valid(*)
1369 integer, intent(in) :: points_id
1370
1371 end subroutine yac_fset_lmask
1372
1373 end interface
1374
1375 !----------------------------------------------------------------------
1379 !----------------------------------------------------------------------
1380
1382
1383 subroutine yac_fdef_imask ( grid_id, &
1384 nbr_points, &
1385 location, &
1386 is_valid, &
1387 mask_id )
1388
1389 integer, intent(in) :: grid_id
1390 integer, intent(in) :: nbr_points
1391 integer, intent(in) :: location
1392 integer, intent(in) :: is_valid(*)
1395 integer, intent(out) :: mask_id
1396
1397 end subroutine yac_fdef_imask
1398
1399 subroutine yac_fdef_lmask ( grid_id, &
1400 nbr_points, &
1401 location, &
1402 is_valid, &
1403 mask_id )
1404
1405 integer, intent(in) :: grid_id
1406 integer, intent(in) :: nbr_points
1407 integer, intent(in) :: location
1408 logical, intent(in) :: is_valid(*)
1411 integer, intent(out) :: mask_id
1412
1413 end subroutine yac_fdef_lmask
1414
1415 end interface
1416
1418
1419 subroutine yac_fdef_imask_named ( grid_id, &
1420 nbr_points, &
1421 location, &
1422 is_valid, &
1423 name, &
1424 mask_id )
1425
1426 integer, intent(in) :: grid_id
1427 integer, intent(in) :: nbr_points
1428 integer, intent(in) :: location
1429 integer, intent(in) :: is_valid(*)
1432 character(len=*), intent(in) :: name
1433 integer, intent(out) :: mask_id
1434
1435 end subroutine yac_fdef_imask_named
1436
1437 subroutine yac_fdef_lmask_named ( grid_id, &
1438 nbr_points, &
1439 location, &
1440 is_valid, &
1441 name, &
1442 mask_id )
1443
1444 integer, intent(in) :: grid_id
1445 integer, intent(in) :: nbr_points
1446 integer, intent(in) :: location
1447 logical, intent(in) :: is_valid(*)
1450 character(len=*), intent(in) :: name
1451 integer, intent(out) :: mask_id
1452
1453 end subroutine yac_fdef_lmask_named
1454
1455 end interface
1456
1457 !----------------------------------------------------------------------
1462 !----------------------------------------------------------------------
1463
1464 interface
1465
1466 subroutine yac_fdef_field ( field_name, &
1467 component_id, &
1468 point_ids, &
1469 num_pointsets, &
1470 collection_size,&
1471 timestep, &
1472 time_unit, &
1473 field_id )
1474
1475 character(len=*), intent (in) :: field_name
1476 integer, intent (in) :: component_id
1477 integer, intent (in) :: point_ids(*)
1478 integer, intent (in) :: num_pointsets
1479 integer, intent (in) :: collection_size
1480 character(len=*), intent (in) :: timestep
1481 integer, intent (in) :: time_unit
1482 integer, intent (out) :: field_id
1483
1484 end subroutine yac_fdef_field
1485
1486 end interface
1487
1488 !----------------------------------------------------------------------
1493 !----------------------------------------------------------------------
1494
1495 interface
1496
1497 subroutine yac_fdef_field_mask ( field_name, &
1498 component_id, &
1499 point_ids, &
1500 mask_ids, &
1501 num_pointsets, &
1502 collection_size,&
1503 timestep, &
1504 time_unit, &
1505 field_id )
1506
1507 character(len=*), intent (in) :: field_name
1508 integer, intent (in) :: component_id
1509 integer, intent (in) :: point_ids(*)
1510 integer, intent (in) :: mask_ids(*)
1511 integer, intent (in) :: num_pointsets
1512 integer, intent (in) :: collection_size
1513 character(len=*), intent (in) :: timestep
1514 integer, intent (in) :: time_unit
1515 integer, intent (out) :: field_id
1516
1517 end subroutine yac_fdef_field_mask
1518
1519 end interface
1520
1521 !----------------------------------------------------------------------
1525 !----------------------------------------------------------------------
1526
1527 interface
1528
1529 subroutine yac_fcheck_field_dimensions( field_id, &
1530 collection_size, &
1531 num_interp_fields, &
1532 interp_field_sizes )
1533
1534 integer, intent (in) :: field_id
1535 integer, intent (in) :: collection_size
1536 integer, intent (in) :: num_interp_fields
1538 integer, intent (in) :: interp_field_sizes(num_interp_fields)
1540
1541 end subroutine yac_fcheck_field_dimensions
1542
1543 end interface
1544
1545 !----------------------------------------------------------------------
1549 !----------------------------------------------------------------------
1550
1551 interface
1552
1553 subroutine yac_fcheck_src_field_buffer_size( field_id, &
1554 collection_size, &
1555 src_field_buffer_size )
1556
1557 integer, intent (in) :: field_id
1558 integer, intent (in) :: collection_size
1559 integer, intent (in) :: src_field_buffer_size
1561
1563
1565 num_src_fields, &
1566 collection_size, &
1567 src_field_buffer_sizes )
1568
1569 integer, intent (in) :: field_id
1570 integer, intent (in) :: num_src_fields
1571 integer, intent (in) :: collection_size
1572 integer, intent (in) :: src_field_buffer_sizes(num_src_fields)
1573
1575
1577
1578 end interface
1579
1580 !----------------------------------------------------------------------
1585 !----------------------------------------------------------------------
1586
1587 interface
1588
1619 field_id, &
1620 frac_mask_fallback_value, &
1621 scaling_factor, &
1622 scaling_summand, &
1623 num_fixed_values, &
1624 fixed_values, &
1625 num_tgt_per_fixed_value, &
1626 tgt_idx_fixed, &
1627 num_wgt_tgt, &
1628 wgt_tgt_idx, &
1629 num_src_per_tgt, &
1630 weights, &
1631 src_field_idx, &
1632 src_idx, &
1633 num_src_fields, &
1634 src_field_buffer_sizes )
1635
1636 integer, intent (in) :: field_id
1637 double precision, intent(out) :: frac_mask_fallback_value
1638 double precision, intent(out) :: scaling_factor
1639 double precision, intent(out) :: scaling_summand
1640 integer, intent(out) :: num_fixed_values
1641 double precision, allocatable, intent(out) :: fixed_values(:)
1642 integer, allocatable, intent(out) :: num_tgt_per_fixed_value(:)
1643 integer, allocatable, intent(out) :: tgt_idx_fixed(:)
1644 integer, intent(out) :: num_wgt_tgt
1645 integer, allocatable, intent(out) :: wgt_tgt_idx(:)
1646 integer, allocatable, intent(out) :: num_src_per_tgt(:)
1647 double precision, allocatable, intent(out) :: weights(:)
1648 integer, allocatable, intent(out) :: src_field_idx(:)
1649 integer, allocatable, intent(out) :: src_idx(:)
1650 integer, intent(out) :: num_src_fields
1651 integer, allocatable, intent(out) :: src_field_buffer_sizes(:)
1653
1682 field_id, &
1683 frac_mask_fallback_value, &
1684 scaling_factor, &
1685 scaling_summand, &
1686 num_fixed_values, &
1687 fixed_values, &
1688 num_tgt_per_fixed_value, &
1689 tgt_idx_fixed, &
1690 src_indptr, &
1691 weights, &
1692 src_field_idx, &
1693 src_idx, &
1694 num_src_fields, &
1695 src_field_buffer_sizes )
1696
1697 integer, intent (in) :: field_id
1698 double precision, intent(out) :: frac_mask_fallback_value
1699 double precision, intent(out) :: scaling_factor
1700 double precision, intent(out) :: scaling_summand
1701 integer, intent(out) :: num_fixed_values
1702 double precision, allocatable, intent(out) :: fixed_values(:)
1703 integer, allocatable, intent(out) :: num_tgt_per_fixed_value(:)
1704 integer, allocatable, intent(out) :: tgt_idx_fixed(:)
1705 integer, allocatable, intent(out) :: src_indptr(:)
1706 double precision, allocatable, intent(out) :: weights(:)
1707 integer, allocatable, intent(out) :: src_field_idx(:)
1708 integer, allocatable, intent(out) :: src_idx(:)
1709 integer, intent(out) :: num_src_fields
1710 integer, allocatable, intent(out) :: src_field_buffer_sizes(:)
1712 end interface
1713
1714 !----------------------------------------------------------------------
1718 !----------------------------------------------------------------------
1719
1720 interface yac_fput
1721
1722 subroutine yac_fput_real ( field_id, &
1723 nbr_hor_points, &
1724 nbr_pointsets, &
1725 collection_size, &
1726 send_field, &
1727 info, &
1728 ierror )
1729
1730 integer, intent (in) :: field_id
1731 integer, intent (in) :: nbr_hor_points
1732 integer, intent (in) :: nbr_pointsets
1733 integer, intent (in) :: collection_size
1734 real, intent (in) :: send_field(nbr_hor_points, &
1735 nbr_pointsets, &
1736 collection_size)
1737
1738 integer, intent (out) :: info
1739 integer, intent (out) :: ierror
1740
1741 end subroutine yac_fput_real
1742
1743 subroutine yac_fput_real_ptr ( field_id, &
1744 nbr_pointsets, &
1745 collection_size, &
1746 send_field, &
1747 info, &
1748 ierror )
1749
1750 import :: yac_real_ptr
1751
1752 integer, intent (in) :: field_id
1753 integer, intent (in) :: nbr_pointsets
1754 integer, intent (in) :: collection_size
1755 type(yac_real_ptr), intent (in) :: send_field(nbr_pointsets, collection_size)
1756
1757 integer, intent (out) :: info
1758 integer, intent (out) :: ierror
1759
1760 end subroutine yac_fput_real_ptr
1761
1762 subroutine yac_fput_single_pointset_real ( field_id, &
1763 nbr_hor_points, &
1764 collection_size, &
1765 send_field, &
1766 info, &
1767 ierror )
1768
1769 integer, intent (in) :: field_id
1770 integer, intent (in) :: nbr_hor_points
1771 integer, intent (in) :: collection_size
1772 real, intent (in) :: send_field(nbr_hor_points, collection_size)
1773
1774 integer, intent (out) :: info
1775 integer, intent (out) :: ierror
1776
1777 end subroutine yac_fput_single_pointset_real
1778
1779 subroutine yac_fput_dble ( field_id, &
1780 nbr_hor_points, &
1781 nbr_pointsets, &
1782 collection_size, &
1783 send_field, &
1784 info, &
1785 ierror )
1786
1787 integer, intent (in) :: field_id
1788 integer, intent (in) :: nbr_hor_points
1789 integer, intent (in) :: nbr_pointsets
1790 integer, intent (in) :: collection_size
1791 double precision, intent (in) :: send_field(nbr_hor_points, nbr_pointsets, collection_size)
1792
1793 integer, intent (out) :: info
1794 integer, intent (out) :: ierror
1795
1796 end subroutine yac_fput_dble
1797
1798 subroutine yac_fput_dble_ptr ( field_id, &
1799 nbr_pointsets, &
1800 collection_size, &
1801 send_field, &
1802 info, &
1803 ierror )
1804
1805 import :: yac_dble_ptr
1806
1807 integer, intent (in) :: field_id
1808 integer, intent (in) :: nbr_pointsets
1809 integer, intent (in) :: collection_size
1810 type(yac_dble_ptr), intent (in) :: send_field(nbr_pointsets, collection_size)
1811
1812 integer, intent (out) :: info
1813 integer, intent (out) :: ierror
1814
1815 end subroutine yac_fput_dble_ptr
1816
1817 subroutine yac_fput_single_pointset_dble ( field_id, &
1818 nbr_hor_points, &
1819 collection_size, &
1820 send_field, &
1821 info, &
1822 ierror )
1823
1824 integer, intent (in) :: field_id
1825 integer, intent (in) :: nbr_hor_points
1826 integer, intent (in) :: collection_size
1827 double precision, intent (in) :: send_field(nbr_hor_points, collection_size)
1828
1829 integer, intent (out) :: info
1830 integer, intent (out) :: ierror
1831
1832 end subroutine yac_fput_single_pointset_dble
1833
1834 subroutine yac_fput_frac_real ( field_id, &
1835 nbr_hor_points, &
1836 nbr_pointsets, &
1837 collection_size, &
1838 send_field, &
1839 send_frac_mask, &
1840 info, &
1841 ierror )
1842
1843 integer, intent (in) :: field_id
1844 integer, intent (in) :: nbr_hor_points
1845 integer, intent (in) :: nbr_pointsets
1846 integer, intent (in) :: collection_size
1847 real, intent (in) :: send_field(nbr_hor_points, &
1848 nbr_pointsets, &
1849 collection_size)
1850
1851 real, intent (in) :: send_frac_mask(nbr_hor_points, &
1852 nbr_pointsets, &
1853 collection_size)
1854
1855 integer, intent (out) :: info
1856 integer, intent (out) :: ierror
1857
1858 end subroutine yac_fput_frac_real
1859
1860 subroutine yac_fput_frac_real_ptr ( field_id, &
1861 nbr_pointsets, &
1862 collection_size, &
1863 send_field, &
1864 send_frac_mask, &
1865 info, &
1866 ierror )
1867
1868 import :: yac_real_ptr
1869
1870 integer, intent (in) :: field_id
1871 integer, intent (in) :: nbr_pointsets
1872 integer, intent (in) :: collection_size
1873 type(yac_real_ptr), intent (in) :: send_field(nbr_pointsets, collection_size)
1874
1875 type(yac_real_ptr), intent (in) :: send_frac_mask(nbr_pointsets, collection_size)
1876
1877 integer, intent (out) :: info
1878 integer, intent (out) :: ierror
1879
1880 end subroutine yac_fput_frac_real_ptr
1881
1882 subroutine yac_fput_frac_single_pointset_real ( field_id, &
1883 nbr_hor_points, &
1884 collection_size, &
1885 send_field, &
1886 send_frac_mask, &
1887 info, &
1888 ierror )
1889
1890 integer, intent (in) :: field_id
1891 integer, intent (in) :: nbr_hor_points
1892 integer, intent (in) :: collection_size
1893 real, intent (in) :: send_field(nbr_hor_points, collection_size)
1894
1895 real, intent (in) :: send_frac_mask(nbr_hor_points, collection_size)
1896
1897 integer, intent (out) :: info
1898 integer, intent (out) :: ierror
1899
1901
1902 subroutine yac_fput_frac_dble ( field_id, &
1903 nbr_hor_points, &
1904 nbr_pointsets, &
1905 collection_size, &
1906 send_field, &
1907 send_frac_mask, &
1908 info, &
1909 ierror )
1910
1911 integer, intent (in) :: field_id
1912 integer, intent (in) :: nbr_hor_points
1913 integer, intent (in) :: nbr_pointsets
1914 integer, intent (in) :: collection_size
1915 double precision, intent (in) :: send_field(nbr_hor_points, nbr_pointsets, collection_size)
1916
1917 double precision, intent (in) :: send_frac_mask(nbr_hor_points, nbr_pointsets, collection_size)
1918
1919 integer, intent (out) :: info
1920 integer, intent (out) :: ierror
1921
1922 end subroutine yac_fput_frac_dble
1923
1924 subroutine yac_fput_frac_dble_ptr ( field_id, &
1925 nbr_pointsets, &
1926 collection_size, &
1927 send_field, &
1928 send_frac_mask, &
1929 info, &
1930 ierror )
1931
1932 import :: yac_dble_ptr
1933
1934 integer, intent (in) :: field_id
1935 integer, intent (in) :: nbr_pointsets
1936 integer, intent (in) :: collection_size
1937 type(yac_dble_ptr), intent (in) :: send_field(nbr_pointsets, collection_size)
1938
1939 type(yac_dble_ptr), intent (in) :: send_frac_mask(nbr_pointsets, collection_size)
1940
1941 integer, intent (out) :: info
1942 integer, intent (out) :: ierror
1943
1944 end subroutine yac_fput_frac_dble_ptr
1945
1946 subroutine yac_fput_frac_single_pointset_dble ( field_id, &
1947 nbr_hor_points, &
1948 collection_size, &
1949 send_field, &
1950 send_frac_mask, &
1951 info, &
1952 ierror )
1953
1954 integer, intent (in) :: field_id
1955 integer, intent (in) :: nbr_hor_points
1956 integer, intent (in) :: collection_size
1957 double precision, intent (in) :: send_field(nbr_hor_points, collection_size)
1958
1959 double precision, intent (in) :: send_frac_mask(nbr_hor_points, collection_size)
1960
1961 integer, intent (out) :: info
1962 integer, intent (out) :: ierror
1963
1965
1966 end interface yac_fput
1967
1968 !----------------------------------------------------------------------
1972 !----------------------------------------------------------------------
1973
1974 interface yac_fget
1975
1976 subroutine yac_fget_real ( field_id, &
1977 nbr_hor_points, &
1978 collection_size, &
1979 recv_field, &
1980 info, &
1981 ierror )
1982
1983 integer, intent (in) :: field_id
1984 integer, intent (in) :: nbr_hor_points
1985 integer, intent (in) :: collection_size
1986 real, intent (inout) :: recv_field(nbr_hor_points, collection_size)
1987 integer, intent (out) :: info
1988 integer, intent (out) :: ierror
1989
1990 end subroutine yac_fget_real
1991
1992 subroutine yac_fget_real_ptr ( field_id, &
1993 collection_size, &
1994 recv_field, &
1995 info, &
1996 ierror )
1997
1998 import :: yac_real_ptr
1999
2000 integer, intent (in) :: field_id
2001 integer, intent (in) :: collection_size
2002 type(yac_real_ptr) :: recv_field(collection_size)
2003 integer, intent (out) :: info
2004 integer, intent (out) :: ierror
2005
2006 end subroutine yac_fget_real_ptr
2007
2008 subroutine yac_fget_dble ( field_id, &
2009 nbr_hor_points, &
2010 collection_size, &
2011 recv_field, &
2012 info, &
2013 ierror )
2014
2015 integer, intent (in) :: field_id
2016 integer, intent (in) :: nbr_hor_points
2017 integer, intent (in) :: collection_size
2018 double precision, intent (inout) :: &
2019 recv_field(nbr_hor_points, collection_size)
2020 integer, intent (out) :: info
2021 integer, intent (out) :: ierror
2022
2023 end subroutine yac_fget_dble
2024
2025 subroutine yac_fget_dble_ptr ( field_id, &
2026 collection_size, &
2027 recv_field, &
2028 info, &
2029 ierror )
2030
2031 import :: yac_dble_ptr
2032
2033 integer, intent (in) :: field_id
2034 integer, intent (in) :: collection_size
2035 type(yac_dble_ptr) :: recv_field(collection_size)
2036 integer, intent (out) :: info
2037 integer, intent (out) :: ierror
2038
2039 end subroutine yac_fget_dble_ptr
2040
2041 end interface yac_fget
2042
2043 !----------------------------------------------------------------------
2047 !----------------------------------------------------------------------
2048
2050
2051 subroutine yac_fget_raw_real ( field_id, &
2052 src_field_buffer_size, &
2053 collection_size, &
2054 src_field_buffer, &
2055 info, &
2056 ierror )
2057
2058 integer, intent (in) :: field_id
2059 integer, intent (in) :: src_field_buffer_size
2061 integer, intent (in) :: collection_size
2062 real, intent (out) :: &
2063 src_field_buffer(src_field_buffer_size, collection_size)
2064
2065 integer, intent (out) :: info
2066 integer, intent (out) :: ierror
2067
2068 end subroutine yac_fget_raw_real
2069
2070 subroutine yac_fget_raw_real_ptr ( field_id, &
2071 num_src_fields, &
2072 collection_size, &
2073 src_field_buffer, &
2074 info, &
2075 ierror )
2076
2077 import :: yac_real_ptr
2078
2079 integer, intent (in) :: field_id
2080 integer, intent (in) :: num_src_fields
2081 integer, intent (in) :: collection_size
2082 type(yac_real_ptr), intent(inout) :: &
2083 src_field_buffer(num_src_fields, collection_size)
2084
2085 integer, intent (out) :: info
2086 integer, intent (out) :: ierror
2087
2088 end subroutine yac_fget_raw_real_ptr
2089
2090 subroutine yac_fget_raw_dble ( field_id, &
2091 src_field_buffer_size, &
2092 collection_size, &
2093 src_field_buffer, &
2094 info, &
2095 ierror )
2096
2097 integer, intent (in) :: field_id
2098 integer, intent (in) :: src_field_buffer_size
2100 integer, intent (in) :: collection_size
2101 double precision, intent (out) :: &
2102 src_field_buffer(src_field_buffer_size, collection_size)
2103
2104 integer, intent (out) :: info
2105 integer, intent (out) :: ierror
2106
2107 end subroutine yac_fget_raw_dble
2108
2109 subroutine yac_fget_raw_dble_ptr ( field_id, &
2110 num_src_fields, &
2111 collection_size, &
2112 src_field_buffer, &
2113 info, &
2114 ierror )
2115
2116 import :: yac_dble_ptr
2117
2118 integer, intent (in) :: field_id
2119 integer, intent (in) :: num_src_fields
2120 integer, intent (in) :: collection_size
2121 type(yac_dble_ptr), intent(inout) :: &
2122 src_field_buffer(num_src_fields, collection_size)
2123
2124 integer, intent (out) :: info
2125 integer, intent (out) :: ierror
2126
2127 end subroutine yac_fget_raw_dble_ptr
2128
2129 subroutine yac_fget_raw_frac_real ( field_id, &
2130 src_field_buffer_size, &
2131 collection_size, &
2132 src_field_buffer, &
2133 src_frac_mask_buffer, &
2134 info, &
2135 ierror )
2136
2137 integer, intent (in) :: field_id
2138 integer, intent (in) :: src_field_buffer_size
2140 integer, intent (in) :: collection_size
2141 real, intent (out) :: &
2142 src_field_buffer(src_field_buffer_size, collection_size)
2143
2144 real, intent (out) :: &
2145 src_frac_mask_buffer(src_field_buffer_size, collection_size)
2146
2147 integer, intent (out) :: info
2148 integer, intent (out) :: ierror
2149
2150 end subroutine yac_fget_raw_frac_real
2151
2152 subroutine yac_fget_raw_frac_real_ptr ( field_id, &
2153 num_src_fields, &
2154 collection_size, &
2155 src_field_buffer, &
2156 src_frac_mask_buffer, &
2157 info, &
2158 ierror )
2159
2160 import :: yac_real_ptr
2161
2162 integer, intent (in) :: field_id
2163 integer, intent (in) :: num_src_fields
2164 integer, intent (in) :: collection_size
2165 type(yac_real_ptr), intent(inout) :: &
2166 src_field_buffer(num_src_fields, collection_size)
2167
2168 type(yac_real_ptr), intent(inout) :: &
2169 src_frac_mask_buffer(num_src_fields, collection_size)
2170
2171 integer, intent (out) :: info
2172 integer, intent (out) :: ierror
2173
2174 end subroutine yac_fget_raw_frac_real_ptr
2175
2176 subroutine yac_fget_raw_frac_dble ( field_id, &
2177 src_field_buffer_size, &
2178 collection_size, &
2179 src_field_buffer, &
2180 src_frac_mask_buffer, &
2181 info, &
2182 ierror )
2183
2184 integer, intent (in) :: field_id
2185 integer, intent (in) :: src_field_buffer_size
2187 integer, intent (in) :: collection_size
2188 double precision, intent (out) :: &
2189 src_field_buffer(src_field_buffer_size, collection_size)
2190
2191 double precision, intent (out) :: &
2192 src_frac_mask_buffer(src_field_buffer_size, collection_size)
2193
2194 integer, intent (out) :: info
2195 integer, intent (out) :: ierror
2196
2197 end subroutine yac_fget_raw_frac_dble
2198
2199 subroutine yac_fget_raw_frac_dble_ptr ( field_id, &
2200 num_src_fields, &
2201 collection_size, &
2202 src_field_buffer, &
2203 src_frac_mask_buffer, &
2204 info, &
2205 ierror )
2206
2207 import :: yac_dble_ptr
2208
2209 integer, intent (in) :: field_id
2210 integer, intent (in) :: num_src_fields
2211 integer, intent (in) :: collection_size
2212 type(yac_dble_ptr), intent(inout) :: &
2213 src_field_buffer(num_src_fields, collection_size)
2214
2215 type(yac_dble_ptr), intent(inout) :: &
2216 src_frac_mask_buffer(num_src_fields, collection_size)
2217
2218 integer, intent (out) :: info
2219 integer, intent (out) :: ierror
2220
2221 end subroutine yac_fget_raw_frac_dble_ptr
2222
2223 end interface yac_fget_raw
2224
2225 !----------------------------------------------------------------------
2229 !----------------------------------------------------------------------
2230
2232
2233 subroutine yac_fget_async_dble_ptr ( field_id, &
2234 collection_size, &
2235 recv_field, &
2236 info, &
2237 ierror )
2238
2239 import :: yac_dble_ptr
2240
2241 integer, intent (in) :: field_id
2242 integer, intent (in) :: collection_size
2243 type(yac_dble_ptr) :: recv_field(collection_size)
2244 integer, intent (out) :: info
2245 integer, intent (out) :: ierror
2246
2247 end subroutine yac_fget_async_dble_ptr
2248
2249 end interface yac_fget_async
2250
2251 !----------------------------------------------------------------------
2256 !----------------------------------------------------------------------
2257
2259
2260 subroutine yac_fget_raw_async_dble_ptr ( field_id, &
2261 num_src_fields, &
2262 collection_size, &
2263 src_field_buffer, &
2264 info, &
2265 ierror )
2266
2267 import :: yac_dble_ptr
2268
2269 integer, intent (in) :: field_id
2270 integer, intent (in) :: num_src_fields
2271 integer, intent (in) :: collection_size
2272 type(yac_dble_ptr), intent(inout) :: &
2273 src_field_buffer(num_src_fields, collection_size)
2274
2275 integer, intent (out) :: info
2276 integer, intent (out) :: ierror
2277
2278 end subroutine yac_fget_raw_async_dble_ptr
2279
2280 subroutine yac_fget_raw_frac_async_dble_ptr ( field_id, &
2281 num_src_fields, &
2282 collection_size, &
2283 src_field_buffer, &
2284 src_frac_mask_buffer, &
2285 info, &
2286 ierror )
2287
2288 import :: yac_dble_ptr
2289
2290 integer, intent (in) :: field_id
2291 integer, intent (in) :: num_src_fields
2292 integer, intent (in) :: collection_size
2293 type(yac_dble_ptr), intent(inout) :: &
2294 src_field_buffer(num_src_fields, collection_size)
2295
2296 type(yac_dble_ptr), intent(inout) :: &
2297 src_frac_mask_buffer(num_src_fields, collection_size)
2298
2299 integer, intent (out) :: info
2300 integer, intent (out) :: ierror
2301
2303
2304 end interface yac_fget_raw_async
2305
2306 !----------------------------------------------------------------------
2310 !----------------------------------------------------------------------
2311
2313
2314 subroutine yac_fexchange_real ( send_field_id, &
2315 recv_field_id, &
2316 send_nbr_hor_points, &
2317 send_nbr_pointsets, &
2318 recv_nbr_hor_points, &
2319 collection_size, &
2320 send_field, &
2321 recv_field, &
2322 send_info, &
2323 recv_info, &
2324 ierror )
2325
2326 integer, intent (in) :: send_field_id
2327 integer, intent (in) :: recv_field_id
2328 integer, intent (in) :: send_nbr_hor_points
2329 integer, intent (in) :: send_nbr_pointsets
2330 integer, intent (in) :: recv_nbr_hor_points
2331 integer, intent (in) :: collection_size
2332 real, intent (in) :: send_field(send_nbr_hor_points, &
2333 send_nbr_pointsets, &
2334 collection_size)
2335
2336 real, intent (inout) :: recv_field(recv_nbr_hor_points, &
2337 collection_size)
2338
2339 integer, intent (out) :: send_info
2340 integer, intent (out) :: recv_info
2341 integer, intent (out) :: ierror
2342
2343 end subroutine yac_fexchange_real
2344
2345 subroutine yac_fexchange_real_ptr ( send_field_id, &
2346 recv_field_id, &
2347 send_nbr_pointsets, &
2348 collection_size, &
2349 send_field, &
2350 recv_field, &
2351 send_info, &
2352 recv_info, &
2353 ierror )
2354
2355 import :: yac_real_ptr
2356
2357 integer, intent (in) :: send_field_id
2358 integer, intent (in) :: recv_field_id
2359 integer, intent (in) :: send_nbr_pointsets
2360 integer, intent (in) :: collection_size
2361 type(yac_real_ptr), intent (in) :: &
2362 send_field(send_nbr_pointsets, &
2363 collection_size)
2364
2365 type(yac_real_ptr) :: recv_field(collection_size)
2366
2367 integer, intent (out) :: send_info
2368 integer, intent (out) :: recv_info
2369 integer, intent (out) :: ierror
2370
2371 end subroutine yac_fexchange_real_ptr
2372
2373 subroutine yac_fexchange_single_pointset_real ( send_field_id, &
2374 recv_field_id, &
2375 send_nbr_hor_points, &
2376 recv_nbr_hor_points, &
2377 collection_size, &
2378 send_field, &
2379 recv_field, &
2380 send_info, &
2381 recv_info, &
2382 ierror )
2383
2384 integer, intent (in) :: send_field_id
2385 integer, intent (in) :: recv_field_id
2386 integer, intent (in) :: send_nbr_hor_points
2387 integer, intent (in) :: recv_nbr_hor_points
2388 integer, intent (in) :: collection_size
2389 real, intent (in) :: send_field(send_nbr_hor_points, &
2390 collection_size)
2391
2392 real, intent (inout) :: recv_field(recv_nbr_hor_points, &
2393 collection_size)
2394
2395 integer, intent (out) :: send_info
2396 integer, intent (out) :: recv_info
2397 integer, intent (out) :: ierror
2398
2400
2401 subroutine yac_fexchange_dble ( send_field_id, &
2402 recv_field_id, &
2403 send_nbr_hor_points, &
2404 send_nbr_pointsets, &
2405 recv_nbr_hor_points, &
2406 collection_size, &
2407 send_field, &
2408 recv_field, &
2409 send_info, &
2410 recv_info, &
2411 ierror )
2412
2413 integer, intent (in) :: send_field_id
2414 integer, intent (in) :: recv_field_id
2415 integer, intent (in) :: send_nbr_hor_points
2416 integer, intent (in) :: send_nbr_pointsets
2417 integer, intent (in) :: recv_nbr_hor_points
2418 integer, intent (in) :: collection_size
2419 double precision, intent (in) :: &
2420 send_field(send_nbr_hor_points, &
2421 send_nbr_pointsets, &
2422 collection_size)
2423
2424 double precision, intent (inout):: &
2425 recv_field(recv_nbr_hor_points, &
2426 collection_size)
2427
2428 integer, intent (out) :: send_info
2429 integer, intent (out) :: recv_info
2430 integer, intent (out) :: ierror
2431
2432 end subroutine yac_fexchange_dble
2433
2434 subroutine yac_fexchange_dble_ptr ( send_field_id, &
2435 recv_field_id, &
2436 send_nbr_pointsets, &
2437 collection_size, &
2438 send_field, &
2439 recv_field, &
2440 send_info, &
2441 recv_info, &
2442 ierror )
2443
2444 import :: yac_dble_ptr
2445
2446 integer, intent (in) :: send_field_id
2447 integer, intent (in) :: recv_field_id
2448 integer, intent (in) :: send_nbr_pointsets
2449 integer, intent (in) :: collection_size
2450 type(yac_dble_ptr), intent (in) :: &
2451 send_field(send_nbr_pointsets, &
2452 collection_size)
2453
2454 type(yac_dble_ptr) :: recv_field(collection_size)
2455
2456 integer, intent (out) :: send_info
2457 integer, intent (out) :: recv_info
2458 integer, intent (out) :: ierror
2459
2460 end subroutine yac_fexchange_dble_ptr
2461
2462 subroutine yac_fexchange_single_pointset_dble ( send_field_id, &
2463 recv_field_id, &
2464 send_nbr_hor_points, &
2465 recv_nbr_hor_points, &
2466 collection_size, &
2467 send_field, &
2468 recv_field, &
2469 send_info, &
2470 recv_info, &
2471 ierror )
2472
2473 integer, intent (in) :: send_field_id
2474 integer, intent (in) :: recv_field_id
2475 integer, intent (in) :: send_nbr_hor_points
2476 integer, intent (in) :: recv_nbr_hor_points
2477 integer, intent (in) :: collection_size
2478 double precision, intent (in) :: &
2479 send_field(send_nbr_hor_points, &
2480 collection_size)
2481
2482 double precision, intent (inout):: &
2483 recv_field(recv_nbr_hor_points, &
2484 collection_size)
2485
2486 integer, intent (out) :: send_info
2487 integer, intent (out) :: recv_info
2488 integer, intent (out) :: ierror
2489
2491
2492 subroutine yac_fexchange_frac_real ( send_field_id, &
2493 recv_field_id, &
2494 send_nbr_hor_points, &
2495 send_nbr_pointsets, &
2496 recv_nbr_hor_points, &
2497 collection_size, &
2498 send_field, &
2499 send_frac_mask, &
2500 recv_field, &
2501 send_info, &
2502 recv_info, &
2503 ierror )
2504
2505 integer, intent (in) :: send_field_id
2506 integer, intent (in) :: recv_field_id
2507 integer, intent (in) :: send_nbr_hor_points
2508 integer, intent (in) :: send_nbr_pointsets
2509 integer, intent (in) :: recv_nbr_hor_points
2510 integer, intent (in) :: collection_size
2511 real, intent (in) :: send_field(send_nbr_hor_points, &
2512 send_nbr_pointsets, &
2513 collection_size)
2514
2515 real, intent (in) :: send_frac_mask(send_nbr_hor_points, &
2516 send_nbr_pointsets, &
2517 collection_size)
2518
2519 real, intent (inout) :: recv_field(recv_nbr_hor_points, &
2520 collection_size)
2521
2522 integer, intent (out) :: send_info
2523 integer, intent (out) :: recv_info
2524 integer, intent (out) :: ierror
2525
2526 end subroutine yac_fexchange_frac_real
2527
2528 subroutine yac_fexchange_frac_real_ptr ( send_field_id, &
2529 recv_field_id, &
2530 send_nbr_pointsets, &
2531 collection_size, &
2532 send_field, &
2533 send_frac_mask, &
2534 recv_field, &
2535 send_info, &
2536 recv_info, &
2537 ierror )
2538
2539 import :: yac_real_ptr
2540
2541 integer, intent (in) :: send_field_id
2542 integer, intent (in) :: recv_field_id
2543 integer, intent (in) :: send_nbr_pointsets
2544 integer, intent (in) :: collection_size
2545 type(yac_real_ptr), intent (in) :: &
2546 send_field(send_nbr_pointsets, &
2547 collection_size)
2548
2549 type(yac_real_ptr), intent (in) :: &
2550 send_frac_mask(send_nbr_pointsets, &
2551 collection_size)
2552
2553 type(yac_real_ptr) :: recv_field(collection_size)
2554
2555 integer, intent (out) :: send_info
2556 integer, intent (out) :: recv_info
2557 integer, intent (out) :: ierror
2558
2559 end subroutine yac_fexchange_frac_real_ptr
2560
2561 subroutine yac_fexchange_frac_single_pointset_real ( send_field_id, &
2562 recv_field_id, &
2563 send_nbr_hor_points, &
2564 recv_nbr_hor_points, &
2565 collection_size, &
2566 send_field, &
2567 send_frac_mask, &
2568 recv_field, &
2569 send_info, &
2570 recv_info, &
2571 ierror )
2572
2573 integer, intent (in) :: send_field_id
2574 integer, intent (in) :: recv_field_id
2575 integer, intent (in) :: send_nbr_hor_points
2576 integer, intent (in) :: recv_nbr_hor_points
2577 integer, intent (in) :: collection_size
2578 real, intent (in) :: send_field(send_nbr_hor_points, &
2579 collection_size)
2580
2581 real, intent (in) :: send_frac_mask(send_nbr_hor_points, &
2582 collection_size)
2583
2584 real, intent (inout) :: recv_field(recv_nbr_hor_points, &
2585 collection_size)
2586
2587 integer, intent (out) :: send_info
2588 integer, intent (out) :: recv_info
2589 integer, intent (out) :: ierror
2590
2592
2593 subroutine yac_fexchange_frac_dble ( send_field_id, &
2594 recv_field_id, &
2595 send_nbr_hor_points, &
2596 send_nbr_pointsets, &
2597 recv_nbr_hor_points, &
2598 collection_size, &
2599 send_field, &
2600 send_frac_mask, &
2601 recv_field, &
2602 send_info, &
2603 recv_info, &
2604 ierror )
2605
2606 integer, intent (in) :: send_field_id
2607 integer, intent (in) :: recv_field_id
2608 integer, intent (in) :: send_nbr_hor_points
2609 integer, intent (in) :: send_nbr_pointsets
2610 integer, intent (in) :: recv_nbr_hor_points
2611 integer, intent (in) :: collection_size
2612 double precision, intent (in) :: &
2613 send_field(send_nbr_hor_points, &
2614 send_nbr_pointsets, &
2615 collection_size)
2616
2617 double precision, intent (in) :: &
2618 send_frac_mask(send_nbr_hor_points, &
2619 send_nbr_pointsets, &
2620 collection_size)
2621
2622 double precision, intent (inout):: &
2623 recv_field(recv_nbr_hor_points, &
2624 collection_size)
2625
2626 integer, intent (out) :: send_info
2627 integer, intent (out) :: recv_info
2628 integer, intent (out) :: ierror
2629
2630 end subroutine yac_fexchange_frac_dble
2631
2632 subroutine yac_fexchange_frac_dble_ptr ( send_field_id, &
2633 recv_field_id, &
2634 send_nbr_pointsets, &
2635 collection_size, &
2636 send_field, &
2637 send_frac_mask, &
2638 recv_field, &
2639 send_info, &
2640 recv_info, &
2641 ierror )
2642
2643 import :: yac_dble_ptr
2644
2645 integer, intent (in) :: send_field_id
2646 integer, intent (in) :: recv_field_id
2647 integer, intent (in) :: send_nbr_pointsets
2648 integer, intent (in) :: collection_size
2649 type(yac_dble_ptr), intent (in) :: &
2650 send_field(send_nbr_pointsets, &
2651 collection_size)
2652
2653 type(yac_dble_ptr), intent (in) :: &
2654 send_frac_mask(send_nbr_pointsets, &
2655 collection_size)
2656
2657 type(yac_dble_ptr) :: recv_field(collection_size)
2658
2659 integer, intent (out) :: send_info
2660 integer, intent (out) :: recv_info
2661 integer, intent (out) :: ierror
2662
2663 end subroutine yac_fexchange_frac_dble_ptr
2664
2665 subroutine yac_fexchange_frac_single_pointset_dble ( send_field_id, &
2666 recv_field_id, &
2667 send_nbr_hor_points, &
2668 recv_nbr_hor_points, &
2669 collection_size, &
2670 send_field, &
2671 send_frac_mask, &
2672 recv_field, &
2673 send_info, &
2674 recv_info, &
2675 ierror )
2676
2677 integer, intent (in) :: send_field_id
2678 integer, intent (in) :: recv_field_id
2679 integer, intent (in) :: send_nbr_hor_points
2680 integer, intent (in) :: recv_nbr_hor_points
2681 integer, intent (in) :: collection_size
2682 double precision, intent (in) :: &
2683 send_field(send_nbr_hor_points, &
2684 collection_size)
2685
2686 double precision, intent (in) :: &
2687 send_frac_mask(send_nbr_hor_points, &
2688 collection_size)
2689
2690 double precision, intent (inout):: &
2691 recv_field(recv_nbr_hor_points, &
2692 collection_size)
2693
2694 integer, intent (out) :: send_info
2695 integer, intent (out) :: recv_info
2696 integer, intent (out) :: ierror
2697
2699
2700 end interface yac_fexchange
2701
2702 !----------------------------------------------------------------------
2707 !----------------------------------------------------------------------
2708
2710
2711 subroutine yac_fexchange_raw_real ( send_field_id, &
2712 recv_field_id, &
2713 send_nbr_hor_points, &
2714 send_nbr_pointsets, &
2715 src_field_buffer_size, &
2716 collection_size, &
2717 send_field, &
2718 src_field_buffer, &
2719 send_info, &
2720 recv_info, &
2721 ierror )
2722
2723 integer, intent (in) :: send_field_id
2724 integer, intent (in) :: recv_field_id
2725 integer, intent (in) :: send_nbr_hor_points
2726 integer, intent (in) :: send_nbr_pointsets
2727 integer, intent (in) :: src_field_buffer_size
2729 integer, intent (in) :: collection_size
2730 real, intent (in) :: &
2731 send_field(send_nbr_hor_points, send_nbr_pointsets, collection_size)
2732
2733 real, intent (out) :: &
2734 src_field_buffer(src_field_buffer_size, collection_size)
2735
2736 integer, intent (out) :: send_info
2737 integer, intent (out) :: recv_info
2738 integer, intent (out) :: ierror
2739
2740 end subroutine yac_fexchange_raw_real
2741
2742 subroutine yac_fexchange_raw_real_ptr ( send_field_id, &
2743 recv_field_id, &
2744 send_nbr_pointsets, &
2745 collection_size, &
2746 send_field, &
2747 src_field_buffer, &
2748 send_info, &
2749 recv_info, &
2750 ierror )
2751
2752 import :: yac_real_ptr
2753
2754 integer, intent (in) :: send_field_id
2755 integer, intent (in) :: recv_field_id
2756 integer, intent (in) :: send_nbr_pointsets
2757 integer, intent (in) :: collection_size
2758 type(yac_real_ptr), intent (in) :: &
2759 send_field(send_nbr_pointsets, collection_size)
2760 type(yac_real_ptr), intent (inout) :: &
2761 src_field_buffer(send_nbr_pointsets, collection_size)
2762 integer, intent (out) :: send_info
2763 integer, intent (out) :: recv_info
2764 integer, intent (out) :: ierror
2765
2766 end subroutine yac_fexchange_raw_real_ptr
2767
2768 subroutine yac_fexchange_raw_single_pointset_real ( send_field_id, &
2769 recv_field_id, &
2770 send_nbr_hor_points, &
2771 src_field_buffer_size, &
2772 collection_size, &
2773 send_field, &
2774 src_field_buffer, &
2775 send_info, &
2776 recv_info, &
2777 ierror )
2778
2779 integer, intent (in) :: send_field_id
2780 integer, intent (in) :: recv_field_id
2781 integer, intent (in) :: send_nbr_hor_points
2782 integer, intent (in) :: src_field_buffer_size
2783 integer, intent (in) :: collection_size
2784 real, intent (in) :: send_field(send_nbr_hor_points, &
2785 collection_size)
2786
2787 real, intent (out) :: &
2788 src_field_buffer(src_field_buffer_size, collection_size)
2789
2790 integer, intent (out) :: send_info
2791 integer, intent (out) :: recv_info
2792 integer, intent (out) :: ierror
2793
2795
2796 subroutine yac_fexchange_raw_dble ( send_field_id, &
2797 recv_field_id, &
2798 send_nbr_hor_points, &
2799 send_nbr_pointsets, &
2800 src_field_buffer_size, &
2801 collection_size, &
2802 send_field, &
2803 src_field_buffer, &
2804 send_info, &
2805 recv_info, &
2806 ierror )
2807
2808 integer, intent (in) :: send_field_id
2809 integer, intent (in) :: recv_field_id
2810 integer, intent (in) :: send_nbr_hor_points
2811 integer, intent (in) :: send_nbr_pointsets
2812 integer, intent (in) :: src_field_buffer_size
2814 integer, intent (in) :: collection_size
2815 double precision, intent (in) :: &
2816 send_field(send_nbr_hor_points, &
2817 send_nbr_pointsets, &
2818 collection_size)
2819
2820 double precision, intent (out) :: &
2821 src_field_buffer(src_field_buffer_size, collection_size)
2822
2823 integer, intent (out) :: send_info
2824 integer, intent (out) :: recv_info
2825 integer, intent (out) :: ierror
2826
2827 end subroutine yac_fexchange_raw_dble
2828
2829 subroutine yac_fexchange_raw_dble_ptr ( send_field_id, &
2830 recv_field_id, &
2831 send_nbr_pointsets, &
2832 collection_size, &
2833 send_field, &
2834 src_field_buffer, &
2835 send_info, &
2836 recv_info, &
2837 ierror )
2838
2839 import :: yac_dble_ptr
2840
2841 integer, intent (in) :: send_field_id
2842 integer, intent (in) :: recv_field_id
2843 integer, intent (in) :: send_nbr_pointsets
2844 integer, intent (in) :: collection_size
2845 type(yac_dble_ptr), intent (in) :: &
2846 send_field(send_nbr_pointsets, collection_size)
2847 type(yac_dble_ptr), intent (inout) :: &
2848 src_field_buffer(send_nbr_pointsets, collection_size)
2849 integer, intent (out) :: send_info
2850 integer, intent (out) :: recv_info
2851 integer, intent (out) :: ierror
2852
2853 end subroutine yac_fexchange_raw_dble_ptr
2854
2855 subroutine yac_fexchange_raw_single_pointset_dble ( send_field_id, &
2856 recv_field_id, &
2857 send_nbr_hor_points, &
2858 src_field_buffer_size, &
2859 collection_size, &
2860 send_field, &
2861 src_field_buffer, &
2862 send_info, &
2863 recv_info, &
2864 ierror )
2865
2866 integer, intent (in) :: send_field_id
2867 integer, intent (in) :: recv_field_id
2868 integer, intent (in) :: send_nbr_hor_points
2869 integer, intent (in) :: src_field_buffer_size
2870 integer, intent (in) :: collection_size
2871 double precision, intent (in) :: &
2872 send_field(send_nbr_hor_points, collection_size)
2873 double precision, intent (out) :: &
2874 src_field_buffer(src_field_buffer_size, collection_size)
2875
2876 integer, intent (out) :: send_info
2877 integer, intent (out) :: recv_info
2878 integer, intent (out) :: ierror
2879
2881
2882 subroutine yac_fexchange_raw_frac_real ( send_field_id, &
2883 recv_field_id, &
2884 send_nbr_hor_points, &
2885 send_nbr_pointsets, &
2886 src_field_buffer_size, &
2887 collection_size, &
2888 send_field, &
2889 send_frac_mask, &
2890 src_field_buffer, &
2891 src_frac_mask_buffer, &
2892 send_info, &
2893 recv_info, &
2894 ierror )
2895
2896 integer, intent (in) :: send_field_id
2897 integer, intent (in) :: recv_field_id
2898 integer, intent (in) :: send_nbr_hor_points
2899 integer, intent (in) :: send_nbr_pointsets
2900 integer, intent (in) :: src_field_buffer_size
2902 integer, intent (in) :: collection_size
2903 real, intent (in) :: &
2904 send_field(send_nbr_hor_points, send_nbr_pointsets, collection_size)
2905
2906 real, intent (in) :: &
2907 send_frac_mask(send_nbr_hor_points, send_nbr_pointsets, collection_size)
2908
2909 real, intent (out) :: &
2910 src_field_buffer(src_field_buffer_size, collection_size)
2911
2912 real, intent (out) :: &
2913 src_frac_mask_buffer(src_field_buffer_size, collection_size)
2914
2915 integer, intent (out) :: send_info
2916 integer, intent (out) :: recv_info
2917 integer, intent (out) :: ierror
2918
2919 end subroutine yac_fexchange_raw_frac_real
2920
2921 subroutine yac_fexchange_raw_frac_real_ptr ( send_field_id, &
2922 recv_field_id, &
2923 send_nbr_pointsets, &
2924 collection_size, &
2925 send_field, &
2926 send_frac_mask, &
2927 src_field_buffer, &
2928 src_frac_mask_buffer, &
2929 send_info, &
2930 recv_info, &
2931 ierror )
2932
2933 import :: yac_real_ptr
2934
2935 integer, intent (in) :: send_field_id
2936 integer, intent (in) :: recv_field_id
2937 integer, intent (in) :: send_nbr_pointsets
2938 integer, intent (in) :: collection_size
2939 type(yac_real_ptr), intent (in) :: &
2940 send_field(send_nbr_pointsets, collection_size)
2941 type(yac_real_ptr), intent (in) :: &
2942 send_frac_mask(send_nbr_pointsets, collection_size)
2943 type(yac_real_ptr), intent(inout) :: &
2944 src_field_buffer(send_nbr_pointsets, collection_size)
2945 type(yac_real_ptr), intent(inout) :: &
2946 src_frac_mask_buffer(send_nbr_pointsets, collection_size)
2947
2948 integer, intent (out) :: send_info
2949 integer, intent (out) :: recv_info
2950 integer, intent (out) :: ierror
2951
2952 end subroutine yac_fexchange_raw_frac_real_ptr
2953
2955 recv_field_id, &
2956 send_nbr_hor_points, &
2957 src_field_buffer_size, &
2958 collection_size, &
2959 send_field, &
2960 send_frac_mask, &
2961 src_field_buffer, &
2962 src_frac_mask_buffer, &
2963 send_info, &
2964 recv_info, &
2965 ierror )
2966
2967 integer, intent (in) :: send_field_id
2968 integer, intent (in) :: recv_field_id
2969 integer, intent (in) :: send_nbr_hor_points
2970 integer, intent (in) :: src_field_buffer_size
2972 integer, intent (in) :: collection_size
2973 real, intent (in) :: &
2974 send_field(send_nbr_hor_points, collection_size)
2975
2976 real, intent (in) :: &
2977 send_frac_mask(send_nbr_hor_points, collection_size)
2978
2979 real, intent (out) :: &
2980 src_field_buffer(src_field_buffer_size, collection_size)
2981
2982 real, intent (out) :: &
2983 src_frac_mask_buffer(src_field_buffer_size, collection_size)
2984
2985 integer, intent (out) :: send_info
2986 integer, intent (out) :: recv_info
2987 integer, intent (out) :: ierror
2988
2990
2991 subroutine yac_fexchange_raw_frac_dble ( send_field_id, &
2992 recv_field_id, &
2993 send_nbr_hor_points, &
2994 send_nbr_pointsets, &
2995 src_field_buffer_size, &
2996 collection_size, &
2997 send_field, &
2998 send_frac_mask, &
2999 src_field_buffer, &
3000 src_frac_mask_buffer, &
3001 send_info, &
3002 recv_info, &
3003 ierror )
3004
3005 integer, intent (in) :: send_field_id
3006 integer, intent (in) :: recv_field_id
3007 integer, intent (in) :: send_nbr_hor_points
3008 integer, intent (in) :: send_nbr_pointsets
3009 integer, intent (in) :: src_field_buffer_size
3011 integer, intent (in) :: collection_size
3012 double precision, intent (in) :: &
3013 send_field(send_nbr_hor_points, send_nbr_pointsets, collection_size)
3014
3015 double precision, intent (in) :: &
3016 send_frac_mask(send_nbr_hor_points, send_nbr_pointsets, collection_size)
3017
3018 double precision, intent (out):: &
3019 src_field_buffer(src_field_buffer_size, collection_size)
3020
3021 double precision, intent (out):: &
3022 src_frac_mask_buffer(src_field_buffer_size, collection_size)
3023
3024 integer, intent (out) :: send_info
3025 integer, intent (out) :: recv_info
3026 integer, intent (out) :: ierror
3027
3028 end subroutine yac_fexchange_raw_frac_dble
3029
3030 subroutine yac_fexchange_raw_frac_dble_ptr ( send_field_id, &
3031 recv_field_id, &
3032 send_nbr_pointsets, &
3033 collection_size, &
3034 send_field, &
3035 send_frac_mask, &
3036 src_field_buffer, &
3037 src_frac_mask_buffer, &
3038 send_info, &
3039 recv_info, &
3040 ierror )
3041
3042 import :: yac_dble_ptr
3043
3044 integer, intent (in) :: send_field_id
3045 integer, intent (in) :: recv_field_id
3046 integer, intent (in) :: send_nbr_pointsets
3047 integer, intent (in) :: collection_size
3048 type(yac_dble_ptr), intent (in) :: &
3049 send_field(send_nbr_pointsets, collection_size)
3050 type(yac_dble_ptr), intent (in) :: &
3051 send_frac_mask(send_nbr_pointsets, collection_size)
3052 type(yac_dble_ptr), intent(inout) :: &
3053 src_field_buffer(send_nbr_pointsets, collection_size)
3054 type(yac_dble_ptr), intent(inout) :: &
3055 src_frac_mask_buffer(send_nbr_pointsets, collection_size)
3056
3057 integer, intent (out) :: send_info
3058 integer, intent (out) :: recv_info
3059 integer, intent (out) :: ierror
3060
3061 end subroutine yac_fexchange_raw_frac_dble_ptr
3062
3064 recv_field_id, &
3065 send_nbr_hor_points, &
3066 src_field_buffer_size, &
3067 collection_size, &
3068 send_field, &
3069 send_frac_mask, &
3070 src_field_buffer, &
3071 src_frac_mask_buffer, &
3072 send_info, &
3073 recv_info, &
3074 ierror )
3075
3076 integer, intent (in) :: send_field_id
3077 integer, intent (in) :: recv_field_id
3078 integer, intent (in) :: send_nbr_hor_points
3079 integer, intent (in) :: src_field_buffer_size
3081 integer, intent (in) :: collection_size
3082 double precision, intent (in) :: &
3083 send_field(send_nbr_hor_points, collection_size)
3084 double precision, intent (in) :: &
3085 send_frac_mask(send_nbr_hor_points, collection_size)
3086 double precision, intent (out) :: &
3087 src_field_buffer(src_field_buffer_size, collection_size)
3088
3089 double precision, intent (out) :: &
3090 src_frac_mask_buffer(src_field_buffer_size, collection_size)
3091
3092 integer, intent (out) :: send_info
3093 integer, intent (out) :: recv_info
3094 integer, intent (out) :: ierror
3095
3097
3098 end interface yac_fexchange_raw
3099
3100 !----------------------------------------------------------------------
3104 !----------------------------------------------------------------------
3105
3106 interface yac_ftest
3107
3108 subroutine yac_ftest_i ( field_id, flag )
3109
3110 integer, intent (in) :: field_id
3111 integer, intent (out) :: flag
3114
3115 end subroutine yac_ftest_i
3116
3117 subroutine yac_ftest_l ( field_id, flag )
3118
3119 integer, intent (in) :: field_id
3120 logical, intent (out) :: flag
3123
3124 end subroutine yac_ftest_l
3125
3126 end interface yac_ftest
3127
3128 !----------------------------------------------------------------------
3132 !----------------------------------------------------------------------
3133
3134 interface yac_fwait
3135
3136 subroutine yac_fwait ( field_id )
3137
3138 integer, intent (in) :: field_id
3139
3140 end subroutine yac_fwait
3141
3142 end interface yac_fwait
3143
3144 !----------------------------------------------------------------------
3148 !----------------------------------------------------------------------
3149
3150 interface
3151
3152 subroutine yac_fget_comp_comm ( comp_id, comp_comm )
3153
3154 integer, intent (in) :: comp_id
3155 integer, intent (out) :: comp_comm
3156
3157 end subroutine yac_fget_comp_comm
3158
3159 end interface
3160
3161 !----------------------------------------------------------------------
3166 !----------------------------------------------------------------------
3167
3169
3170 subroutine yac_fget_comps_comm ( comp_names, &
3171 num_comps, &
3172 comps_comm)
3173
3174 use, intrinsic :: iso_c_binding, only : c_char
3175 import :: yac_max_charlen
3176
3177 integer, intent(in) :: num_comps
3178 character(kind=c_char, len=YAC_MAX_CHARLEN), intent(in) :: &
3179 comp_names(num_comps)
3180 integer, intent (out) :: comps_comm
3181
3182 end subroutine yac_fget_comps_comm
3183
3184 subroutine yac_fget_comps_comm_instance ( yac_instance_id, &
3185 comp_names, &
3186 num_comps, &
3187 comps_comm)
3188
3189 use, intrinsic :: iso_c_binding, only : c_char
3190 import :: yac_max_charlen
3191
3192 integer, intent(in) :: yac_instance_id
3193 integer, intent(in) :: num_comps
3194 character(kind=c_char, len=YAC_MAX_CHARLEN), intent(in) :: &
3195 comp_names(num_comps)
3196 integer, intent (out) :: comps_comm
3197
3198 end subroutine yac_fget_comps_comm_instance
3199
3200 end interface yac_fget_comps_comm
3201
3202
3203 !----------------------------------------------------------------------
3207 !----------------------------------------------------------------------
3208
3210
3211 subroutine yac_fsync_def ( )
3212
3213 end subroutine yac_fsync_def
3214
3215 subroutine yac_fsync_def_instance ( yac_instance_id )
3216
3217 integer, intent(in) :: yac_instance_id
3218
3219 end subroutine yac_fsync_def_instance
3220
3221 end interface yac_fsync_def
3222
3223 !----------------------------------------------------------------------
3227 !----------------------------------------------------------------------
3228
3229 interface
3230
3231 subroutine yac_fget_interp_stack_config(interp_stack_config_id)
3232 integer, intent(out) :: interp_stack_config_id
3233 end subroutine yac_fget_interp_stack_config
3234
3235 subroutine yac_ffree_interp_stack_config(interp_stack_config_id)
3236 integer, intent(in) :: interp_stack_config_id
3237 end subroutine yac_ffree_interp_stack_config
3238
3239 subroutine yac_fadd_interp_stack_config_average(interp_stack_config_id, &
3240 reduction_type, partial_coverage)
3241 integer, intent(in) :: interp_stack_config_id
3242 integer, intent(in) :: reduction_type
3243 integer, intent(in) :: partial_coverage
3245
3246 subroutine yac_fadd_interp_stack_config_ncc(interp_stack_config_id, &
3247 weight_type, partial_coverage)
3248 integer, intent(in) :: interp_stack_config_id
3249 integer, intent(in) :: weight_type
3250 integer, intent(in) :: partial_coverage
3252
3253 subroutine yac_fadd_interp_stack_config_nnn(interp_stack_config_id, &
3254 type, n, max_search_distance, scale)
3255 integer, intent(in) :: interp_stack_config_id
3256 integer, intent(in) :: type
3257 integer, intent(in) :: n
3258 double precision, intent(in) :: max_search_distance
3259 double precision, intent(in) :: scale
3261
3262 subroutine yac_fadd_interp_stack_config_rbf(interp_stack_config_id, &
3263 n, max_search_distance, scale)
3264 integer, intent(in) :: interp_stack_config_id
3265 integer, optional, intent(in) :: n
3266 double precision, optional, intent(in) :: max_search_distance
3267 double precision, optional, intent(in) :: scale
3269
3270 subroutine yac_fadd_interp_stack_config_conservative(interp_stack_config_id, &
3271 order, enforced_conserv, partial_coverage, normalization)
3272 integer, intent(in) :: interp_stack_config_id
3273 integer, intent(in) :: order
3274 integer, intent(in) :: enforced_conserv
3275 integer, intent(in) :: partial_coverage
3276 integer, intent(in) :: normalization
3278
3279 subroutine yac_fadd_interp_stack_config_spmap(interp_stack_config_id, &
3280 spread_distance, max_search_distance, weight_type, scale_type, &
3281 src_sphere_radius, src_filename, src_varname, src_min_global_id, &
3282 tgt_sphere_radius, tgt_filename, tgt_varname, tgt_min_global_id)
3283 integer, intent(in) :: interp_stack_config_id
3284 double precision, intent(in) :: spread_distance
3285 double precision, intent(in) :: max_search_distance
3286 integer, intent(in) :: weight_type
3287 integer, intent(in) :: scale_type
3288 double precision, intent(in) :: src_sphere_radius
3289 character (len=*), intent(in) :: src_filename
3290 character (len=*), intent(in) :: src_varname
3291 integer, intent(in) :: src_min_global_id
3292 double precision, intent(in) :: tgt_sphere_radius
3293 character (len=*), intent(in) :: tgt_filename
3294 character (len=*), intent(in) :: tgt_varname
3295 integer, intent(in) :: tgt_min_global_id
3297
3298 subroutine yac_fadd_interp_stack_config_hcsbb(interp_stack_config_id)
3299 integer, intent(in) :: interp_stack_config_id
3301
3303 interp_stack_config_id, filename, on_missing_file, on_success)
3304 integer, intent(in) :: interp_stack_config_id
3305 character (len=*), intent(in) :: filename
3306 integer, optional, intent(in) :: on_missing_file
3307 integer, optional, intent(in) :: on_success
3309
3310 subroutine yac_fadd_interp_stack_config_fixed(interp_stack_config_id, &
3311 val)
3312 integer, intent(in) :: interp_stack_config_id
3313 double precision, intent(in) :: val
3315
3316 subroutine yac_fadd_interp_stack_config_creep(interp_stack_config_id, &
3317 creep_distance)
3318 integer, intent(in) :: interp_stack_config_id
3319 integer, intent(in) :: creep_distance
3321
3323 interp_stack_config, interp_stack_config_id)
3324 character ( len=* ), intent(in) :: interp_stack_config
3325 integer, intent(out) :: interp_stack_config_id
3327
3329 interp_stack_config, interp_stack_config_id)
3330 character ( len=* ), intent(in) :: interp_stack_config
3331 integer, intent(out) :: interp_stack_config_id
3333
3334 end interface
3335
3336 !----------------------------------------------------------------------
3340 !----------------------------------------------------------------------
3341
3343
3344 subroutine yac_fdef_couple ( &
3345 src_comp_name, src_grid_name, src_field_name, &
3346 tgt_comp_name, tgt_grid_name, tgt_field_name, &
3347 coupling_timestep, time_unit, time_reduction, &
3348 interp_stack_config_id, src_lag, tgt_lag, &
3349 weight_file, weight_file_on_existing, &
3350 mapping_side, scale_factor, scale_summand, &
3351 src_mask_names, tgt_mask_name, &
3352 yaxt_exchanger_name, use_raw_exchange )
3353
3354 import :: yac_string
3355 character ( len=* ), intent(in) :: src_comp_name
3356 character ( len=* ), intent(in) :: src_grid_name
3357 character ( len=* ), intent(in) :: src_field_name
3358 character ( len=* ), intent(in) :: tgt_comp_name
3359 character ( len=* ), intent(in) :: tgt_grid_name
3360 character ( len=* ), intent(in) :: tgt_field_name
3361 character ( len=* ), intent(in) :: coupling_timestep
3362 integer, intent(in) :: time_unit
3363 integer, intent(in) :: time_reduction
3364 integer, intent(in) :: interp_stack_config_id
3365 integer, intent(in), optional :: src_lag
3366 integer, intent(in), optional :: tgt_lag
3367 character ( len=* ), intent(in), optional :: weight_file
3368 integer, intent(in), optional :: weight_file_on_existing
3369 integer, intent(in), optional :: mapping_side
3370 double precision, intent(in), optional :: scale_factor
3371 double precision, intent(in), optional :: scale_summand
3372 type(yac_string), intent(in), optional :: src_mask_names(:)
3373 character ( len=* ), intent(in), optional :: tgt_mask_name
3374 character ( len=* ), intent(in), optional :: yaxt_exchanger_name
3375 logical, intent(in), optional :: use_raw_exchange
3376 end subroutine yac_fdef_couple
3377
3378 subroutine yac_fdef_couple_instance ( instance_id, &
3379 src_comp_name, src_grid_name, src_field_name, &
3380 tgt_comp_name, tgt_grid_name, tgt_field_name, &
3381 coupling_timestep, time_unit, time_reduction, &
3382 interp_stack_config_id, src_lag, tgt_lag, &
3383 weight_file, weight_file_on_existing, &
3384 mapping_side, scale_factor, scale_summand, &
3385 src_mask_names, tgt_mask_name, &
3386 yaxt_exchanger_name, use_raw_exchange )
3387
3388 import :: yac_string
3389 integer, intent(in) :: instance_id
3390 character ( len=* ), intent(in) :: src_comp_name
3391 character ( len=* ), intent(in) :: src_grid_name
3392 character ( len=* ), intent(in) :: src_field_name
3393 character ( len=* ), intent(in) :: tgt_comp_name
3394 character ( len=* ), intent(in) :: tgt_grid_name
3395 character ( len=* ), intent(in) :: tgt_field_name
3396 character ( len=* ), intent(in) :: coupling_timestep
3397 integer, intent(in) :: time_unit
3398 integer, intent(in) :: time_reduction
3399 integer, intent(in) :: interp_stack_config_id
3400 integer, intent(in), optional :: src_lag
3401 integer, intent(in), optional :: tgt_lag
3402 character ( len=* ), intent(in), optional :: weight_file
3403 integer, intent(in), optional :: weight_file_on_existing
3404 integer, intent(in), optional :: mapping_side
3405 double precision, intent(in), optional :: scale_factor
3406 double precision, intent(in), optional :: scale_summand
3407 type(yac_string), intent(in), optional :: src_mask_names(:)
3408 character ( len=* ), intent(in), optional :: tgt_mask_name
3409 character ( len=* ), intent(in), optional :: yaxt_exchanger_name
3410 logical, intent(in), optional :: use_raw_exchange
3411 end subroutine yac_fdef_couple_instance
3412
3413 end interface yac_fdef_couple
3414
3415 !----------------------------------------------------------------------
3419 !----------------------------------------------------------------------
3420
3421 interface yac_fenddef
3422
3423 subroutine yac_fenddef ( )
3424
3425 end subroutine yac_fenddef
3426
3427 subroutine yac_fenddef_instance ( yac_instance_id )
3428
3429 integer, intent(in) :: yac_instance_id
3430
3431 end subroutine yac_fenddef_instance
3432
3433 subroutine yac_fenddef_and_emit_config(emit_flags, config)
3434
3435 integer, intent (in) :: emit_flags
3436 character (len=:), ALLOCATABLE :: config
3437
3438 end subroutine yac_fenddef_and_emit_config
3439
3441 yac_instance_id, emit_flags, config)
3442
3443 integer, intent (in) :: yac_instance_id
3444 integer, intent (in) :: emit_flags
3445 character (len=:), ALLOCATABLE :: config
3446
3448
3449 end interface yac_fenddef
3450
3451 !----------------------------------------------------------------------
3455 !----------------------------------------------------------------------
3456
3458
3459 function yac_fget_grid_size( location, grid_id ) result (grid_size)
3460
3461 integer, intent(in) :: location
3462 integer, intent(in) :: grid_id
3463 integer :: grid_size
3464
3465 end function yac_fget_grid_size
3466
3467 end interface yac_fget_grid_size
3468
3469 ! ---------------------------------------------------------------------
3470
3472
3474 grid_id, nbr_cells, cell_areas )
3475
3476 integer, intent(in) :: grid_id
3477 integer, intent(in) :: nbr_cells
3478 real, intent(out) :: cell_areas(nbr_cells)
3480
3482 grid_id, nbr_cells, cell_areas )
3483
3484 integer, intent(in) :: grid_id
3485 integer, intent(in) :: nbr_cells
3486 double precision, intent(out) :: cell_areas(nbr_cells)
3488
3489 end interface
3490
3491 ! ---------------------------------------------------------------------
3492
3494
3495 function yac_fget_points_size( point_id ) result (points_size)
3496
3497 integer, intent(in) :: point_id
3498 integer :: points_size
3499
3500 end function yac_fget_points_size
3501
3502 end interface yac_fget_points_size
3503
3504 ! ---------------------------------------------------------------------
3505
3507
3508 function yac_fget_field_id ( comp_name, grid_name, field_name ) result(field_id)
3509
3510 character(len=*), intent(in) :: comp_name
3511 character(len=*), intent(in) :: grid_name
3512 character(len=*), intent(in) :: field_name
3513 integer :: field_id
3514
3515 end function yac_fget_field_id
3516
3517 function yac_fget_field_id_instance ( yac_id, &
3518 comp_name, &
3519 grid_name, &
3520 field_name ) result(field_id)
3521
3522 integer, intent(in) :: yac_id
3523 character(len=*), intent(in) :: comp_name
3524 character(len=*), intent(in) :: grid_name
3525 character(len=*), intent(in) :: field_name
3526 integer :: field_id
3527
3528 end function yac_fget_field_id_instance
3529 end interface yac_fget_field_id
3530
3531 ! ---------------------------------------------------------------------
3532
3534
3536 result( comp_name )
3537
3538 integer, intent (in) :: field_id
3539 character(len=:), ALLOCATABLE :: comp_name
3540
3542
3543 end interface yac_fget_component_name
3544
3546
3547 function yac_fget_grid_name_from_field_id ( field_id ) &
3548 result( grid_name )
3549
3550 integer, intent (in) :: field_id
3551 character(len=:), ALLOCATABLE :: grid_name
3552
3554
3555 end interface yac_fget_grid_name
3556
3558
3560 result( field_name )
3561
3562 integer, intent (in) :: field_id
3563 character(len=:), ALLOCATABLE :: field_name
3564
3566
3567 end interface yac_fget_field_name
3568
3569 ! ---------------------------------------------------------------------
3570
3572
3573 function yac_fget_comp_names( ) result(comp_names)
3574 import :: yac_string
3575 type(yac_string), allocatable :: comp_names(:)
3576 end function yac_fget_comp_names
3577
3578 function yac_fget_comp_names_instance( yac_instance_id) &
3579 result( comp_names )
3580 import :: yac_string
3581 integer, intent(in) :: yac_instance_id
3582 type(yac_string), allocatable :: comp_names(:)
3583 end function yac_fget_comp_names_instance
3584
3585 end interface yac_fget_comp_names
3586
3588
3589 function yac_fget_grid_names( ) result( grid_names )
3590 import :: yac_string
3591 type(yac_string), allocatable :: grid_names(:)
3592 end function yac_fget_grid_names
3593
3594 function yac_fget_grid_names_instance( yac_instance_id ) &
3595 result( grid_names )
3596 import :: yac_string
3597 integer, intent(in) :: yac_instance_id
3598 type(yac_string), allocatable :: grid_names(:)
3599 end function yac_fget_grid_names_instance
3600
3601 end interface yac_fget_grid_names
3602
3604
3605 function yac_fget_comp_grid_names( comp_name ) result( grid_names )
3606 import :: yac_string
3607 character(len=*), intent(in) :: comp_name
3608 type(yac_string), allocatable :: grid_names(:)
3609 end function yac_fget_comp_grid_names
3610
3611 function yac_fget_comp_grid_names_instance( yac_instance_id, comp_name ) &
3612 result( grid_names )
3613 import :: yac_string
3614 integer, intent(in) :: yac_instance_id
3615 character(len=*), intent(in) :: comp_name
3616 type(yac_string), allocatable :: grid_names(:)
3618
3619 end interface yac_fget_comp_grid_names
3620
3622
3623 function yac_fget_field_names( comp_name, grid_name ) result ( field_names )
3624 import :: yac_string
3625 character(len=*), intent(in) :: comp_name
3626 character(len=*), intent(in) :: grid_name
3627 type(yac_string), allocatable :: field_names(:)
3628 end function yac_fget_field_names
3629
3630 function yac_fget_field_names_instance( yac_instance_id, &
3631 comp_name, &
3632 grid_name ) &
3633 result( field_names )
3634 import :: yac_string
3635 integer, intent(in) :: yac_instance_id
3636 character(len=*), intent(in) :: comp_name
3637 character(len=*), intent(in) :: grid_name
3638 type(yac_string), allocatable :: field_names(:)
3640
3641 end interface yac_fget_field_names
3642
3643 ! ---------------------------------------------------------------------
3644
3646
3647 function yac_fget_role_from_field_id ( field_id )
3648
3649 integer, intent (in) :: field_id
3651
3652 end function yac_fget_role_from_field_id
3653
3654 function yac_fget_field_role ( comp_name, grid_name, field_name )
3655 character(len=*), intent(in) :: comp_name
3656 character(len=*), intent(in) :: grid_name
3657 character(len=*), intent(in) :: field_name
3658 integer :: yac_fget_field_role
3659 end function yac_fget_field_role
3660
3661 function yac_fget_field_role_instance ( yac_instance_id, &
3662 comp_name, &
3663 grid_name, &
3664 field_name )
3665 integer, intent(in) :: yac_instance_id
3666 character(len=*), intent(in) :: comp_name
3667 character(len=*), intent(in) :: grid_name
3668 character(len=*), intent(in) :: field_name
3670 end function yac_fget_field_role_instance
3671
3672 end interface yac_fget_field_role
3673
3674 ! ---------------------------------------------------------------------
3675
3677
3678 subroutine yac_fget_field_source ( tgt_comp_name, &
3679 tgt_grid_name, &
3680 tgt_field_name, &
3681 src_comp_name, &
3682 src_grid_name, &
3683 src_field_name )
3684 character(len=*), intent(in) :: tgt_comp_name
3685 character(len=*), intent(in) :: tgt_grid_name
3686 character(len=*), intent(in) :: tgt_field_name
3687 character(len=:), ALLOCATABLE :: src_comp_name
3688 character(len=:), ALLOCATABLE :: src_grid_name
3689 character(len=:), ALLOCATABLE :: src_field_name
3690 end subroutine yac_fget_field_source
3691
3692 subroutine yac_fget_field_source_instance ( yac_instance_id, &
3693 tgt_comp_name, &
3694 tgt_grid_name, &
3695 tgt_field_name, &
3696 src_comp_name, &
3697 src_grid_name, &
3698 src_field_name )
3699 integer, intent(in) :: yac_instance_id
3700 character(len=*), intent(in) :: tgt_comp_name
3701 character(len=*), intent(in) :: tgt_grid_name
3702 character(len=*), intent(in) :: tgt_field_name
3703 character(len=:), ALLOCATABLE :: src_comp_name
3704 character(len=:), ALLOCATABLE :: src_grid_name
3705 character(len=:), ALLOCATABLE :: src_field_name
3706 end subroutine yac_fget_field_source_instance
3707
3708 end interface yac_fget_field_source
3709
3710 ! ---------------------------------------------------------------------
3711
3713
3714
3715 function yac_fget_timestep_from_field_id ( field_id ) result( timestep )
3716 integer, intent (in) :: field_id
3717 character(len=:), ALLOCATABLE :: timestep
3718
3720
3721 function yac_fget_field_timestep ( comp_name, grid_name, field_name ) &
3722 result( timestep )
3723
3724 character(len=*), intent(in) :: comp_name
3725 character(len=*), intent(in) :: grid_name
3726 character(len=*), intent(in) :: field_name
3727 character(len=:), ALLOCATABLE :: timestep
3728
3729 end function yac_fget_field_timestep
3730
3731 function yac_fget_field_timestep_instance ( yac_instance_id, &
3732 comp_name, &
3733 grid_name, &
3734 field_name ) &
3735 result( timestep )
3736
3737 integer, intent(in) :: yac_instance_id
3738 character(len=*), intent(in) :: comp_name
3739 character(len=*), intent(in) :: grid_name
3740 character(len=*), intent(in) :: field_name
3741 character(len=:), ALLOCATABLE :: timestep
3742
3744
3745 end interface yac_fget_field_timestep
3746
3747 ! ---------------------------------------------------------------------
3748
3750
3751 function yac_fget_field_datetime ( field_id ) &
3752 result( datetime )
3753 integer, intent(in) :: field_id
3754 character(len=:), allocatable :: datetime
3755 end function yac_fget_field_datetime
3756
3757 end interface yac_fget_field_datetime
3758
3759 ! ---------------------------------------------------------------------
3760
3762
3764 comp_name, grid_name, field_name, frac_mask_fallback_value )
3765 character(len=*), intent(in) :: comp_name
3766 character(len=*), intent(in) :: grid_name
3767 character(len=*), intent(in) :: field_name
3768 double precision, intent(in) :: frac_mask_fallback_value
3769 end subroutine yac_fenable_field_frac_mask
3770
3772 yac_instance_id, comp_name, grid_name, field_name, &
3773 frac_mask_fallback_value )
3774 integer, intent(in) :: yac_instance_id
3775 character(len=*), intent(in) :: comp_name
3776 character(len=*), intent(in) :: grid_name
3777 character(len=*), intent(in) :: field_name
3778 double precision, intent(in) :: frac_mask_fallback_value
3780
3781 end interface yac_fenable_field_frac_mask
3782
3783 ! ---------------------------------------------------------------------
3784
3786
3788 comp_name, grid_name, field_name )
3789 character(len=*), intent(in) :: comp_name
3790 character(len=*), intent(in) :: grid_name
3791 character(len=*), intent(in) :: field_name
3792 double precision :: yac_fget_field_frac_mask_fallback_value
3794
3796 yac_instance_id, comp_name, grid_name, field_name )
3797 integer, intent(in) :: yac_instance_id
3798 character(len=*), intent(in) :: comp_name
3799 character(len=*), intent(in) :: grid_name
3800 character(len=*), intent(in) :: field_name
3803
3805
3806 ! ---------------------------------------------------------------------
3807
3809
3811 result(collection_size)
3812
3813 integer, intent (in) :: field_id
3814 integer :: collection_size
3815
3817
3818 function yac_fget_field_collection_size ( comp_name, grid_name, field_name )
3819 character(len=*), intent(in) :: comp_name
3820 character(len=*), intent(in) :: grid_name
3821 character(len=*), intent(in) :: field_name
3824
3825 function yac_fget_field_collection_size_instance ( yac_instance_id, &
3826 comp_name, &
3827 grid_name, &
3828 field_name )
3829 integer, intent(in) :: yac_instance_id
3830 character(len=*), intent(in) :: comp_name
3831 character(len=*), intent(in) :: grid_name
3832 character(len=*), intent(in) :: field_name
3835
3836 end interface yac_fget_field_collection_size
3837
3838 ! ---------------------------------------------------------------------
3839
3841 subroutine yac_fdef_component_metadata(comp_name, metadata)
3842 character(len=*), intent(in) :: comp_name
3843 character(len=*), intent(in) :: metadata
3844 end subroutine yac_fdef_component_metadata
3845
3846 subroutine yac_fdef_component_metadata_instance(yac_instance_id, comp_name, &
3847 metadata)
3848 integer, intent(in) :: yac_instance_id
3849 character(len=*), intent(in) :: comp_name
3850 character(len=*), intent(in) :: metadata
3852 end interface yac_fdef_component_metadata
3853
3855 function yac_fcomponent_has_metadata(comp_name) result( has_metadata )
3856 character(len=*), intent(in) :: comp_name
3857 logical :: has_metadata
3858 end function yac_fcomponent_has_metadata
3859
3860 function yac_fcomponent_has_metadata_instance(yac_instance_id, comp_name) &
3861 result( has_metadata )
3862 integer, intent(in) :: yac_instance_id
3863 character(len=*), intent(in) :: comp_name
3864 logical :: has_metadata
3866 end interface yac_fcomponent_has_metadata
3867
3869 function yac_fget_component_metadata(comp_name) result( metadata )
3870 character(len=*), intent(in) :: comp_name
3871 character(len=:), allocatable :: metadata
3872 end function yac_fget_component_metadata
3873
3874 function yac_fget_component_metadata_instance(yac_instance_id, comp_name) &
3875 result( metadata )
3876 integer, intent(in) :: yac_instance_id
3877 character(len=*), intent(in) :: comp_name
3878 character(len=:), allocatable :: metadata
3880 end interface yac_fget_component_metadata
3881
3883 subroutine yac_fdef_grid_metadata(grid_name, metadata)
3884 character(len=*), intent(in) :: grid_name
3885 character(len=*), intent(in) :: metadata
3886 end subroutine yac_fdef_grid_metadata
3887
3888 subroutine yac_fdef_grid_metadata_instance(yac_instance_id, grid_name, metadata)
3889 integer, intent(in) :: yac_instance_id
3890 character(len=*), intent(in) :: grid_name
3891 character(len=*), intent(in) :: metadata
3892 end subroutine yac_fdef_grid_metadata_instance
3893 end interface yac_fdef_grid_metadata
3894
3896 function yac_fgrid_has_metadata(grid_name) result( has_metadata )
3897 character(len=*), intent(in) :: grid_name
3898 logical :: has_metadata
3899 end function yac_fgrid_has_metadata
3900
3901 function yac_fgrid_has_metadata_instance(yac_instance_id, grid_name) &
3902 result( has_metadata )
3903 integer, intent(in) :: yac_instance_id
3904 character(len=*), intent(in) :: grid_name
3905 logical :: has_metadata
3907 end interface yac_fgrid_has_metadata
3908
3910 function yac_fget_grid_metadata(grid_name) result( metadata )
3911 character(len=*), intent(in) :: grid_name
3912 character(len=:), allocatable :: metadata
3913 end function yac_fget_grid_metadata
3914
3915 function yac_fget_grid_metadata_instance(yac_instance_id, grid_name) &
3916 result( metadata )
3917 integer, intent(in) :: yac_instance_id
3918 character(len=*), intent(in) :: grid_name
3919 character(len=:), allocatable :: metadata
3921 end interface yac_fget_grid_metadata
3922
3924 subroutine yac_fdef_field_metadata(comp_name, grid_name, field_name, metadata)
3925 character(len=*), intent(in) :: comp_name
3926 character(len=*), intent(in) :: grid_name
3927 character(len=*), intent(in) :: field_name
3928 character(len=*), intent(in) :: metadata
3929 end subroutine yac_fdef_field_metadata
3930
3931 subroutine yac_fdef_field_metadata_instance(yac_instance_id, comp_name, &
3932 grid_name, field_name, metadata)
3933 integer, intent(in) :: yac_instance_id
3934 character(len=*), intent(in) :: comp_name
3935 character(len=*), intent(in) :: grid_name
3936 character(len=*), intent(in) :: field_name
3937 character(len=*), intent(in) :: metadata
3939 end interface yac_fdef_field_metadata
3940
3942 function yac_ffield_has_metadata(comp_name, grid_name, field_name) &
3943 result( has_metadata )
3944 character(len=*), intent(in) :: comp_name
3945 character(len=*), intent(in) :: grid_name
3946 character(len=*), intent(in) :: field_name
3947 logical :: has_metadata
3948 end function yac_ffield_has_metadata
3949
3951 yac_instance_id, comp_name, grid_name, field_name) &
3952 result( has_metadata )
3953 integer, intent(in) :: yac_instance_id
3954 character(len=*), intent(in) :: comp_name
3955 character(len=*), intent(in) :: grid_name
3956 character(len=*), intent(in) :: field_name
3957 logical :: has_metadata
3959 end interface yac_ffield_has_metadata
3960
3962 function yac_fget_field_metadata(comp_name, grid_name, field_name) &
3963 result( metadata )
3964 character(len=*), intent(in) :: comp_name
3965 character(len=*), intent(in) :: grid_name
3966 character(len=*), intent(in) :: field_name
3967 character(len=:), allocatable :: metadata
3968 end function yac_fget_field_metadata
3969
3970 function yac_fget_field_metadata_instance(yac_instance_id, comp_name, &
3971 grid_name, field_name) &
3972 result( metadata )
3973 integer, intent(in) :: yac_instance_id
3974 character(len=*), intent(in) :: comp_name
3975 character(len=*), intent(in) :: grid_name
3976 character(len=*), intent(in) :: field_name
3977 character(len=:), allocatable :: metadata
3979 end interface yac_fget_field_metadata
3980
3981 ! ---------------------------------------------------------------------
3982
3983 interface
3984
3985 subroutine yac_fget_action ( field_id, action )
3986
3987 integer, intent (in) :: field_id
3988 integer, intent (out) :: action
3994
3995 end subroutine yac_fget_action
3996
3997 end interface
3998
3999 ! ---------------------------------------------------------------------
4000
4001 interface
4002
4003 subroutine yac_fupdate ( field_id )
4004
4005 integer, intent (in) :: field_id
4006
4007 end subroutine yac_fupdate
4008
4009 end interface
4010
4011! --- ISO_C interface -------------------------------------------------
4012
4013 interface
4014
4015 subroutine yac_abort_message ( text, file, line ) &
4016 bind( c, name='yac_abort_message' )
4017
4018 use, intrinsic :: iso_c_binding, only : c_int, c_char
4019
4020 character ( kind=c_char ), dimension(*) :: text
4021 character ( kind=c_char ), dimension(*) :: file
4022 integer ( kind=c_int ), value :: line
4023
4024 end subroutine yac_abort_message
4025
4026 end interface
4027
4028end module yac
Fortran interface for checking the dimensions of a field.
Fortran interface for checking the source field buffer sizes.
Fortran interface for the coupler cleanup before restart.
Fortran interface for definition of a couple.
Fortran interface for the definition of time parameters.
Fortran interface for the definition of coupling fields using explicit masks.
Fortran interface for the definition of coupling fields using default masks.
Fortran interface for the definition of grids.
Fortran interface for the definition of masks.
Fortran interface for the definition of points.
Fortran interface for invoking the end of the definition phase.
Fortran interface for exchanging coupling fields using raw data exchange.
Fortran interface for exchanging coupling fields.
Fortran interface for the coupler termination.
Fortran interface for asynchronous receiving coupling fields.
Fortran interface for getting back a local MPI communicator.
Fortran interface for getting back a MPI communicator for communication between components.
Fortran interface for getting default instance id.
Fortran interface for invoking query functions.
Fortran interfaces for the definition of an interpolation stack.
Fortran interface for asynchronous receiving coupling fields using raw exchange.
Get interpolation information for raw data exchange.
Fortran interface for getting interpolation information for raw data exchange.
Fortran interface for receiving coupling fields using raw exchange.
Fortran interface for getting the start- and end datetime.
Fortran interface for getting the yac version.
Fortran interface for receiving coupling fields.
Fortran interface for the coupler initialisation.
Fortran2C interface for yac collective routines.
Fortran interface for the component definition.
Fortran interface for sending coupling fields.
Fortran interface for the reading of configuration files.
Fortran interface for the writing of configuration files.
Fortran interface for the setting of a grid core masks.
Fortran interface for the setting of grid global ids.
Fortran interface for the writing of grid files.
Fortran interface for the setting of default pointset masks.
Fortran interface for invoking the end of the definition phase.
Fortran interface for testing fields for active communicaitons.
Fortran interface for testing fields for active communicaitons.
@ yac_location_edge
@ yac_location_cell
@ yac_location_corner
@ yac_time_unit_millisecond
@ yac_time_unit_minute
@ yac_time_unit_hour
@ yac_time_unit_year
@ yac_time_unit_second
@ yac_time_unit_iso_format
@ yac_time_unit_month
@ yac_time_unit_day
double precision, parameter yac_frac_mask_no_value
@ yac_action_coupling
data exchange
@ yac_action_out_of_bound
put/get is outside of the valid range
@ yac_action_reduction
data reduction, but data exchange
@ yac_action_put_for_restart
last valid put
@ yac_action_none
no data exchanges
@ yac_action_get_for_restart
last valid get
@ yac_nnn_gauss
@ yac_nnn_dist
@ yac_nnn_rbf
@ yac_nnn_avg
@ yac_ncc_avg
@ yac_ncc_dist
@ yac_file_success_cont
continue weight computation with following interpolation methods
@ yac_file_success_stop
prevents following interpolation method from computating further weights
@ yac_wgt_on_existing_keep
@ yac_wgt_on_existing_overwrite
@ yac_wgt_on_existing_error
@ yac_proleptic_gregorian
@ yac_year_of_365_days
@ yac_calendar_not_set
@ yac_year_of_360_days
integer, parameter yac_max_charlen
Constants.
integer yac_yaml_emitter_default_f
Flag paramters for emitting of coupling configurations.
@ yac_conserv_fracarea
@ yac_conserv_destarea
@ yac_avg_bary
@ yac_avg_arithmetic
@ yac_avg_dist
@ yac_spmap_avg
@ yac_spmap_dist
@ yac_config_output_format_json
@ yac_config_output_format_yaml
@ yac_config_output_sync_loc_sync_def
@ yac_config_output_sync_loc_enddef
@ yac_config_output_sync_loc_def_comp
@ yac_reduction_time_average
@ yac_reduction_time_accumulate
@ yac_reduction_time_minimum
@ yac_reduction_time_none
@ yac_reduction_time_maximum
integer, parameter yac_mpi_fint_kind
@ yac_exchange_type_target
@ yac_exchange_type_source
@ yac_exchange_type_none
@ yac_file_missing_error
abort on missing file
@ yac_file_missing_cont
continue on missing file
integer yac_yaml_emitter_json_f
@ yac_spmap_srcarea
@ yac_spmap_invtgtarea
@ yac_spmap_fracarea
@ yac_spmap_none
pointer data types
subroutine yac_fcompute_grid_cell_areas_real(grid_id, nbr_cells, cell_areas)
subroutine yac_fset_grid_output_file_instance(yac_instance_id, gridname, filename)
subroutine yac_fexchange_raw_frac_dble(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, src_field_buffer_size, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fexchange_raw_frac_real_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fsync_def_instance(yac_instance_id)
subroutine yac_fcompute_grid_cell_areas_dble(grid_id, nbr_cells, cell_areas)
subroutine yac_fset_config_output_file_instance(yac_instance_id, filename, fileformat, sync_location, include_definitions)
character(len=:) function, allocatable yac_fget_field_metadata_instance(yac_instance_id, comp_name, grid_name, field_name)
subroutine yac_fexchange_frac_real_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
integer function yac_fget_field_role_instance(yac_instance_id, comp_name, grid_name, field_name)
subroutine yac_fexchange_raw_single_pointset_real(send_field_id, recv_field_id, send_nbr_hor_points, src_field_buffer_size, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
logical function yac_fgrid_has_metadata_instance(yac_instance_id, grid_name)
subroutine yac_fget_raw_frac_real(field_id, src_field_buffer_size, collection_size, src_field_buffer, src_frac_mask_buffer, info, ierror)
subroutine yac_fenddef_and_emit_config_instance(yac_instance_id, emit_flags, config)
subroutine yac_fexchange_dble(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, recv_nbr_hor_points, collection_size, send_field, recv_field, send_info, recv_info, ierror)
subroutine yac_fget_raw_dble(field_id, src_field_buffer_size, collection_size, src_field_buffer, info, ierror)
subroutine yac_fdef_comp_dummy_instance(yac_instance_id)
subroutine yac_fexchange_raw_real(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, src_field_buffer_size, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
character(len=:) function, allocatable yac_fget_component_name_from_field_id(field_id)
subroutine yac_fget_comps_comm_instance(yac_instance_id, comp_names, num_comps, comps_comm)
subroutine yac_fput_frac_single_pointset_real(field_id, nbr_hor_points, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fdef_grid_unstruct_dble(grid_name, nbr_vertices, nbr_cells, nbr_vertices_per_cell_in, x_vertices, y_vertices, cell_to_vertex_in, grid_id, use_ll_edges)
Definition of a uniform unstructured grid (all cells have the number of vertices)
subroutine yac_fput_single_pointset_dble(field_id, nbr_hor_points, collection_size, send_field, info, ierror)
subroutine yac_fdef_lmask_named(grid_id, nbr_points, location, is_valid, name, mask_id)
type(yac_string) function, dimension(:), allocatable yac_fget_comp_names_instance(yac_instance_id)
subroutine yac_fget_real_ptr(field_id, collection_size, recv_field, info, ierror)
subroutine yac_fget_raw_real(field_id, src_field_buffer_size, collection_size, src_field_buffer, info, ierror)
subroutine yac_fget_raw_dble_ptr(field_id, num_src_fields, collection_size, src_field_buffer, info, ierror)
subroutine yac_fdef_couple_instance(instance_id, src_comp_name, src_grid_name, src_field_name, tgt_comp_name, tgt_grid_name, tgt_field_name, coupling_timestep, time_unit, time_reduction, interp_stack_config_id, src_lag, tgt_lag, weight_file, weight_file_on_existing, mapping_side, scale_factor, scale_summand, src_mask_names, tgt_mask_name, yaxt_exchanger_name, use_raw_exchange)
character(len=:) function, allocatable yac_fget_grid_name_from_field_id(field_id)
subroutine yac_fexchange_raw_dble_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
subroutine yac_fexchange_raw_frac_real(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, src_field_buffer_size, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fexchange_frac_dble_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
subroutine yac_finit_instance(yac_instance_id)
character(len=:) function, allocatable yac_fget_field_timestep_instance(yac_instance_id, comp_name, grid_name, field_name)
subroutine yac_fdef_datetime_instance(yac_instance_id, start_datetime, end_datetime)
subroutine yac_fget_raw_frac_real_ptr(field_id, num_src_fields, collection_size, src_field_buffer, src_frac_mask_buffer, info, ierror)
subroutine yac_fput_frac_dble(field_id, nbr_hor_points, nbr_pointsets, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fput_frac_single_pointset_dble(field_id, nbr_hor_points, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fexchange_frac_real(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, recv_nbr_hor_points, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
subroutine yac_fexchange_frac_dble(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, recv_nbr_hor_points, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
subroutine yac_fdef_grid_metadata_instance(yac_instance_id, grid_name, metadata)
subroutine yac_fput_dble_ptr(field_id, nbr_pointsets, collection_size, send_field, info, ierror)
subroutine yac_fexchange_single_pointset_dble(send_field_id, recv_field_id, send_nbr_hor_points, recv_nbr_hor_points, collection_size, send_field, recv_field, send_info, recv_info, ierror)
subroutine yac_fdef_grid_unstruct_real(grid_name, nbr_vertices, nbr_cells, nbr_vertices_per_cell_in, x_vertices_real, y_vertices_real, cell_to_vertex_in, grid_id, use_ll_edges)
Definition of a uniform unstructured grid (all cells have the number of vertices)
subroutine yac_fdef_grid_reg2d_rot_real(grid_name, nbr_vertices, cyclic, x_vertices_real, y_vertices_real, x_north_pole_real, y_north_pole_real, grid_id)
Definition of a 2d regular rotated grid.
subroutine yac_fcleanup_instance(yac_instance_id)
subroutine yac_fget_dble(field_id, nbr_hor_points, collection_size, recv_field, info, ierror)
integer function yac_fget_field_id_instance(yac_id, comp_name, grid_name, field_name)
subroutine yac_fdef_points_curve2d_dble(grid_id, nbr_points, location, x_points, y_points, point_id)
subroutine yac_fexchange_raw_frac_single_pointset_dble(send_field_id, recv_field_id, send_nbr_hor_points, src_field_buffer_size, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fexchange_raw_frac_dble_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fdef_grid_curve2d_dble(grid_name, nbr_vertices, cyclic, x_vertices, y_vertices, grid_id)
Definition of a 2d curvilinear grid.
type(yac_string) function, dimension(:), allocatable yac_fget_comp_grid_names_instance(yac_instance_id, comp_name)
subroutine yac_fget_raw_frac_async_dble_ptr(field_id, num_src_fields, collection_size, src_field_buffer, src_frac_mask_buffer, info, ierror)
subroutine yac_finit_comm_dummy(world_comm)
type(yac_string) function, dimension(:), allocatable yac_fget_grid_names_instance(yac_instance_id)
subroutine yac_fdef_imask(grid_id, nbr_points, location, is_valid, mask_id)
subroutine yac_fdef_points_reg2d_real(grid_id, nbr_points, location, x_points_real, y_points_real, point_id)
character(len=:) function, allocatable yac_fget_timestep_from_field_id(field_id)
character(len=:) function, allocatable yac_fget_field_name_from_field_id(field_id)
subroutine yac_fput_dble(field_id, nbr_hor_points, nbr_pointsets, collection_size, send_field, info, ierror)
subroutine yac_fput_single_pointset_real(field_id, nbr_hor_points, collection_size, send_field, info, ierror)
subroutine yac_finit_comm_instance(comm, yac_instance_id)
logical function yac_fcomponent_has_metadata_instance(yac_instance_id, comp_name)
subroutine yac_fset_lmask(is_valid, points_id)
subroutine yac_fdef_grid_unstruct_edge_real(grid_name, nbr_vertices, nbr_cells, nbr_edges, nbr_edges_per_cell_in, x_vertices_real, y_vertices_real, cell_to_edge_in, edge_to_vertex_in, grid_id, use_ll_edges)
Definition of a uniform unstructured grid (all cells have the number of vertices) with explicit edge ...
subroutine yac_fdef_imask_named(grid_id, nbr_points, location, is_valid, name, mask_id)
subroutine yac_fdef_grid_curve2d_real(grid_name, nbr_vertices, cyclic, x_vertices_real, y_vertices_real, grid_id)
Definition of a 2d curvilinear grid.
type(yac_string) function, dimension(:), allocatable yac_fget_field_names_instance(yac_instance_id, comp_name, grid_name)
subroutine yac_fdef_comp_instance(yac_instance_id, comp_name, comp_id)
subroutine yac_fget_raw_async_dble_ptr(field_id, num_src_fields, collection_size, src_field_buffer, info, ierror)
subroutine yac_fdef_grid_nonuniform_real(grid_name, nbr_vertices, nbr_cells, nbr_connections, nbr_vertices_per_cell, x_vertices_real, y_vertices_real, cell_to_vertex_in, grid_id, use_ll_edges)
Definition of a non-uniform unstructured grid (cells have varying numbers of vertices)
subroutine yac_fdef_points_reg2d_rot_dble(grid_id, nbr_points, location, x_points, y_points, x_north_pole, y_north_pole, point_id)
subroutine yac_fexchange_frac_single_pointset_dble(send_field_id, recv_field_id, send_nbr_hor_points, recv_nbr_hor_points, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
subroutine yac_ftest_l(field_id, flag)
subroutine yac_fexchange_raw_single_pointset_dble(send_field_id, recv_field_id, send_nbr_hor_points, src_field_buffer_size, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
double precision function yac_fget_field_frac_mask_fallback_value_instance(yac_instance_id, comp_name, grid_name, field_name)
logical function yac_ffield_has_metadata_instance(yac_instance_id, comp_name, grid_name, field_name)
subroutine yac_fput_real_ptr(field_id, nbr_pointsets, collection_size, send_field, info, ierror)
subroutine yac_ffinalize_instance(yac_instance_id)
integer function yac_fget_collection_size_from_field_id(field_id)
subroutine yac_fset_core_imask(is_core, location, grid_id)
subroutine yac_fdef_grid_reg2d_rot_dble(grid_name, nbr_vertices, cyclic, x_vertices, y_vertices, x_north_pole, y_north_pole, grid_id)
Definition of a 2d regular rotated grid.
integer function yac_fget_role_from_field_id(field_id)
subroutine yac_fenable_field_frac_mask_instance(yac_instance_id, comp_name, grid_name, field_name, frac_mask_fallback_value)
subroutine yac_fdef_points_reg2d_rot_real(grid_id, nbr_points, location, x_points_real, y_points_real, x_north_pole_real, y_north_pole_real, point_id)
subroutine yac_fget_raw_frac_dble_ptr(field_id, num_src_fields, collection_size, src_field_buffer, src_frac_mask_buffer, info, ierror)
subroutine yac_fget_async_dble_ptr(field_id, collection_size, recv_field, info, ierror)
character(len=:) function, allocatable yac_fget_end_datetime_instance(yac_instance_id)
subroutine yac_fget_raw_real_ptr(field_id, num_src_fields, collection_size, src_field_buffer, info, ierror)
subroutine yac_fexchange_raw_dble(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, src_field_buffer_size, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
subroutine yac_fread_config_yaml_instance(yac_instance_id, yaml_filename)
subroutine yac_fput_frac_dble_ptr(field_id, nbr_pointsets, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fexchange_real(send_field_id, recv_field_id, send_nbr_hor_points, send_nbr_pointsets, recv_nbr_hor_points, collection_size, send_field, recv_field, send_info, recv_info, ierror)
subroutine yac_fput_frac_real(field_id, nbr_hor_points, nbr_pointsets, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fexchange_real_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, recv_field, send_info, recv_info, ierror)
character(len=:) function, allocatable yac_fget_grid_metadata_instance(yac_instance_id, grid_name)
subroutine yac_ftest_i(field_id, flag)
subroutine yac_fset_core_lmask(is_core, location, grid_id)
subroutine yac_fdef_points_curve2d_real(grid_id, nbr_points, location, x_points_real, y_points_real, point_id)
subroutine yac_fget_real(field_id, nbr_hor_points, collection_size, recv_field, info, ierror)
subroutine yac_fdef_grid_unstruct_edge_dble(grid_name, nbr_vertices, nbr_cells, nbr_edges, nbr_edges_per_cell_in, x_vertices, y_vertices, cell_to_edge_in, edge_to_vertex_in, grid_id, use_ll_edges)
Definition of a uniform unstructured grid (all cells have the number of vertices) with explicit edge ...
subroutine yac_fget_field_source_instance(yac_instance_id, tgt_comp_name, tgt_grid_name, tgt_field_name, src_comp_name, src_grid_name, src_field_name)
subroutine yac_fenddef_and_emit_config(emit_flags, config)
integer function yac_fget_field_collection_size_instance(yac_instance_id, comp_name, grid_name, field_name)
subroutine yac_fexchange_dble_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, recv_field, send_info, recv_info, ierror)
subroutine yac_fput_real(field_id, nbr_hor_points, nbr_pointsets, collection_size, send_field, info, ierror)
subroutine yac_fdef_comps_instance(yac_instance_id, comp_names, num_comps, comp_ids)
subroutine yac_fdef_grid_nonuniform_edge_dble(grid_name, nbr_vertices, nbr_cells, nbr_edges, nbr_connections, nbr_edges_per_cell, x_vertices, y_vertices, cell_to_edge_in, edge_to_vertex_in, grid_id, use_ll_edges)
Definition of a non-uniform unstructured grid (cells have varying numbers of vertices) with explicit ...
subroutine yac_fget_dble_ptr(field_id, collection_size, recv_field, info, ierror)
subroutine yac_fenddef_instance(yac_instance_id)
subroutine yac_fdef_grid_cloud_real(grid_name, nbr_points, x_points_real, y_points_real, grid_id)
Definition of a grid consisting of a cloud of points.
character(len=:) function, allocatable yac_fget_component_metadata_instance(yac_instance_id, comp_name)
subroutine yac_fput_frac_real_ptr(field_id, nbr_pointsets, collection_size, send_field, send_frac_mask, info, ierror)
subroutine yac_fdef_lmask(grid_id, nbr_points, location, is_valid, mask_id)
subroutine yac_fdef_component_metadata_instance(yac_instance_id, comp_name, metadata)
subroutine yac_fset_imask(is_valid, points_id)
subroutine yac_fpredef_comp_instance(yac_instance_id, comp_name, comp_id)
subroutine yac_fexchange_single_pointset_real(send_field_id, recv_field_id, send_nbr_hor_points, recv_nbr_hor_points, collection_size, send_field, recv_field, send_info, recv_info, ierror)
subroutine yac_fdef_grid_nonuniform_edge_real(grid_name, nbr_vertices, nbr_cells, nbr_edges, nbr_connections, nbr_edges_per_cell, x_vertices_real, y_vertices_real, cell_to_edge, edge_to_vertex, grid_id, use_ll_edges)
Definition of a non-uniform unstructured grid (cells have varying numbers of vertices) with explicit ...
subroutine yac_fexchange_frac_single_pointset_real(send_field_id, recv_field_id, send_nbr_hor_points, recv_nbr_hor_points, collection_size, send_field, send_frac_mask, recv_field, send_info, recv_info, ierror)
subroutine yac_fdef_grid_nonuniform_dble(grid_name, nbr_vertices, nbr_cells, nbr_connections, nbr_vertices_per_cell, x_vertices, y_vertices, cell_to_vertex_in, grid_id, use_ll_edges)
Definition of a non-uniform unstructured grid (cells have varying numbers of vertices)
subroutine yac_fdef_field_metadata_instance(yac_instance_id, comp_name, grid_name, field_name, metadata)
subroutine yac_fexchange_raw_frac_single_pointset_real(send_field_id, recv_field_id, send_nbr_hor_points, src_field_buffer_size, collection_size, send_field, send_frac_mask, src_field_buffer, src_frac_mask_buffer, send_info, recv_info, ierror)
subroutine yac_fread_config_json_instance(yac_instance_id, json_filename)
subroutine yac_fexchange_raw_real_ptr(send_field_id, recv_field_id, send_nbr_pointsets, collection_size, send_field, src_field_buffer, send_info, recv_info, ierror)
subroutine yac_fdef_grid_reg2d_dble(grid_name, nbr_vertices, cyclic, x_vertices, y_vertices, grid_id)
Definition of a 2d regular grid.
subroutine yac_fdef_grid_cloud_dble(grid_name, nbr_points, x_points, y_points, grid_id)
Definition of a grid consisting of a cloud of points.
subroutine yac_fdef_points_unstruct_real(grid_id, nbr_points, location, x_points_real, y_points_real, point_id)
subroutine yac_fdef_points_reg2d_dble(grid_id, nbr_points, location, x_points, y_points, point_id)
subroutine yac_fdef_grid_reg2d_real(grid_name, nbr_vertices, cyclic, x_vertices_real, y_vertices_real, grid_id)
Definition of a 2d regular grid.
subroutine yac_fget_raw_frac_dble(field_id, src_field_buffer_size, collection_size, src_field_buffer, src_frac_mask_buffer, info, ierror)
subroutine yac_fdef_points_unstruct_dble(grid_id, nbr_points, location, x_points, y_points, point_id)
character(len=:) function, allocatable yac_fget_start_datetime_instance(yac_instance_id)