-
Notifications
You must be signed in to change notification settings - Fork 10
/
Copy pathbmi.f90
564 lines (495 loc) · 20.8 KB
/
bmi.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
! The Basic Model Interface (BMI) Fortran specification.
!
! This language specification is derived from the Scientific
! Interface Definition Language (SIDL) file bmi.sidl located at
! https://github.com/csdms/bmi.
module bmif_2_0
implicit none
integer, parameter :: BMI_MAX_COMPONENT_NAME = 2048
integer, parameter :: BMI_MAX_VAR_NAME = 2048
integer, parameter :: BMI_MAX_TYPE_NAME = 2048
integer, parameter :: BMI_MAX_UNITS_NAME = 2048
integer, parameter :: BMI_FAILURE = 1
integer, parameter :: BMI_SUCCESS = 0
type, abstract :: bmi
contains
! Initialize, run, finalize (IRF)
procedure(bmif_initialize), deferred :: initialize
procedure(bmif_update), deferred :: update
procedure(bmif_update_until), deferred :: update_until
procedure(bmif_finalize), deferred :: finalize
! Exchange items
procedure(bmif_get_component_name), deferred :: get_component_name
procedure(bmif_get_input_item_count), deferred :: get_input_item_count
procedure(bmif_get_output_item_count), deferred :: get_output_item_count
procedure(bmif_get_input_var_names), deferred :: get_input_var_names
procedure(bmif_get_output_var_names), deferred :: get_output_var_names
! Variable information
procedure(bmif_get_var_grid), deferred :: get_var_grid
procedure(bmif_get_var_type), deferred :: get_var_type
procedure(bmif_get_var_units), deferred :: get_var_units
procedure(bmif_get_var_itemsize), deferred :: get_var_itemsize
procedure(bmif_get_var_nbytes), deferred :: get_var_nbytes
procedure(bmif_get_var_location), deferred :: get_var_location
! Time information
procedure(bmif_get_current_time), deferred :: get_current_time
procedure(bmif_get_start_time), deferred :: get_start_time
procedure(bmif_get_end_time), deferred :: get_end_time
procedure(bmif_get_time_units), deferred :: get_time_units
procedure(bmif_get_time_step), deferred :: get_time_step
! Getters, by type
procedure(bmif_get_value_int), deferred :: get_value_int
procedure(bmif_get_value_float), deferred :: get_value_float
procedure(bmif_get_value_double), deferred :: get_value_double
procedure(bmif_get_value_ptr_int), deferred :: get_value_ptr_int
procedure(bmif_get_value_ptr_float), deferred :: get_value_ptr_float
procedure(bmif_get_value_ptr_double), deferred :: get_value_ptr_double
procedure(bmif_get_value_at_indices_int), deferred :: &
get_value_at_indices_int
procedure(bmif_get_value_at_indices_float), deferred :: &
get_value_at_indices_float
procedure(bmif_get_value_at_indices_double), deferred :: &
get_value_at_indices_double
! Setters, by type
procedure(bmif_set_value_int), deferred :: set_value_int
procedure(bmif_set_value_float), deferred :: set_value_float
procedure(bmif_set_value_double), deferred :: set_value_double
procedure(bmif_set_value_at_indices_int), deferred :: &
set_value_at_indices_int
procedure(bmif_set_value_at_indices_float), deferred :: &
set_value_at_indices_float
procedure(bmif_set_value_at_indices_double), deferred :: &
set_value_at_indices_double
! Grid information
procedure(bmif_get_grid_rank), deferred :: get_grid_rank
procedure(bmif_get_grid_size), deferred :: get_grid_size
procedure(bmif_get_grid_type), deferred :: get_grid_type
! Uniform rectilinear
procedure(bmif_get_grid_shape), deferred :: get_grid_shape
procedure(bmif_get_grid_spacing), deferred :: get_grid_spacing
procedure(bmif_get_grid_origin), deferred :: get_grid_origin
! Non-uniform rectilinear, curvilinear
procedure(bmif_get_grid_x), deferred :: get_grid_x
procedure(bmif_get_grid_y), deferred :: get_grid_y
procedure(bmif_get_grid_z), deferred :: get_grid_z
! Unstructured
procedure(bmif_get_grid_node_count), deferred :: get_grid_node_count
procedure(bmif_get_grid_edge_count), deferred :: get_grid_edge_count
procedure(bmif_get_grid_face_count), deferred :: get_grid_face_count
procedure(bmif_get_grid_edge_nodes), deferred :: get_grid_edge_nodes
procedure(bmif_get_grid_face_edges), deferred :: get_grid_face_edges
procedure(bmif_get_grid_face_nodes), deferred :: get_grid_face_nodes
procedure(bmif_get_grid_nodes_per_face), deferred :: &
get_grid_nodes_per_face
end type bmi
abstract interface
! Perform startup tasks for the model.
function bmif_initialize(this, config_file) result(bmi_status)
import :: bmi
class(bmi), intent(out) :: this
character(len=*), intent(in) :: config_file
integer :: bmi_status
end function bmif_initialize
! Advance the model one time step.
function bmif_update(this) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
integer :: bmi_status
end function bmif_update
! Advance the model until the given time.
function bmif_update_until(this, time) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
double precision, intent(in) :: time
integer :: bmi_status
end function bmif_update_until
! Perform teardown tasks for the model.
function bmif_finalize(this) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
integer :: bmi_status
end function bmif_finalize
! Get the name of the model.
function bmif_get_component_name(this, name) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), pointer, intent(out) :: name
integer :: bmi_status
end function bmif_get_component_name
! Count a model's input variables.
function bmif_get_input_item_count(this, count) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(out) :: count
integer :: bmi_status
end function bmif_get_input_item_count
! Count a model's output variables.
function bmif_get_output_item_count(this, count) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(out) :: count
integer :: bmi_status
end function bmif_get_output_item_count
! List a model's input variables.
function bmif_get_input_var_names(this, names) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), pointer, intent(out) :: names(:)
integer :: bmi_status
end function bmif_get_input_var_names
! List a model's output variables.
function bmif_get_output_var_names(this, names) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), pointer, intent(out) :: names(:)
integer :: bmi_status
end function bmif_get_output_var_names
! Get the grid identifier for the given variable.
function bmif_get_var_grid(this, name, grid) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, intent(out) :: grid
integer :: bmi_status
end function bmif_get_var_grid
! Get the data type of the given variable as a string.
function bmif_get_var_type(this, name, type) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
character(len=*), intent(out) :: type
integer :: bmi_status
end function bmif_get_var_type
! Get the units of the given variable.
function bmif_get_var_units(this, name, units) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
character(len=*), intent(out) :: units
integer :: bmi_status
end function bmif_get_var_units
! Get memory use per array element, in bytes.
function bmif_get_var_itemsize(this, name, size) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, intent(out) :: size
integer :: bmi_status
end function bmif_get_var_itemsize
! Get size of the given variable, in bytes.
function bmif_get_var_nbytes(this, name, nbytes) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, intent(out) :: nbytes
integer :: bmi_status
end function bmif_get_var_nbytes
! Describe where a variable is located: node, edge, or face.
function bmif_get_var_location(this, name, location) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
character(len=*), intent(out) :: location
integer :: bmi_status
end function bmif_get_var_location
! Current time of the model.
function bmif_get_current_time(this, time) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
double precision, intent(out) :: time
integer :: bmi_status
end function bmif_get_current_time
! Start time of the model.
function bmif_get_start_time(this, time) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
double precision, intent(out) :: time
integer :: bmi_status
end function bmif_get_start_time
! End time of the model.
function bmif_get_end_time(this, time) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
double precision, intent(out) :: time
integer :: bmi_status
end function bmif_get_end_time
! Time units of the model.
function bmif_get_time_units(this, units) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(out) :: units
integer :: bmi_status
end function bmif_get_time_units
! Time step of the model.
function bmif_get_time_step(this, time_step) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
double precision, intent(out) :: time_step
integer :: bmi_status
end function bmif_get_time_step
! Get a copy of values (flattened!) of the given integer variable.
function bmif_get_value_int(this, name, dest) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, intent(inout) :: dest(:)
integer :: bmi_status
end function bmif_get_value_int
! Get a copy of values (flattened!) of the given real variable.
function bmif_get_value_float(this, name, dest) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
real, intent(inout) :: dest(:)
integer :: bmi_status
end function bmif_get_value_float
! Get a copy of values (flattened!) of the given double variable.
function bmif_get_value_double(this, name, dest) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
double precision, intent(inout) :: dest(:)
integer :: bmi_status
end function bmif_get_value_double
! Get a reference to the given integer variable.
function bmif_get_value_ptr_int(this, name, dest_ptr) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, pointer, intent(inout) :: dest_ptr(:)
integer :: bmi_status
end function bmif_get_value_ptr_int
! Get a reference to the given real variable.
function bmif_get_value_ptr_float(this, name, dest_ptr) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
real, pointer, intent(inout) :: dest_ptr(:)
integer :: bmi_status
end function bmif_get_value_ptr_float
! Get a reference to the given double variable.
function bmif_get_value_ptr_double(this, name, dest_ptr) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
double precision, pointer, intent(inout) :: dest_ptr(:)
integer :: bmi_status
end function bmif_get_value_ptr_double
! Get integer values at particular (one-dimensional) indices.
function bmif_get_value_at_indices_int(this, name, dest, inds) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
integer, intent(inout) :: dest(:)
integer, intent(in) :: inds(:)
integer :: bmi_status
end function bmif_get_value_at_indices_int
! Get real values at particular (one-dimensional) indices.
function bmif_get_value_at_indices_float(this, name, dest, inds) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
real, intent(inout) :: dest(:)
integer, intent(in) :: inds(:)
integer :: bmi_status
end function bmif_get_value_at_indices_float
! Get double values at particular (one-dimensional) indices.
function bmif_get_value_at_indices_double(this, name, dest, inds) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
character(len=*), intent(in) :: name
double precision, intent(inout) :: dest(:)
integer, intent(in) :: inds(:)
integer :: bmi_status
end function bmif_get_value_at_indices_double
! Set new values for an integer model variable.
function bmif_set_value_int(this, name, src) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
integer, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_int
! Set new values for a real model variable.
function bmif_set_value_float(this, name, src) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
real, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_float
! Set new values for a double model variable.
function bmif_set_value_double(this, name, src) result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
double precision, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_double
! Set integer values at particular (one-dimensional) indices.
function bmif_set_value_at_indices_int(this, name, inds, src) &
result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
integer, intent(in) :: inds(:)
integer, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_at_indices_int
! Set real values at particular (one-dimensional) indices.
function bmif_set_value_at_indices_float(this, name, inds, src) &
result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
integer, intent(in) :: inds(:)
real, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_at_indices_float
! Set double values at particular (one-dimensional) indices.
function bmif_set_value_at_indices_double(this, name, inds, src) &
result(bmi_status)
import :: bmi
class(bmi), intent(inout) :: this
character(len=*), intent(in) :: name
integer, intent(in) :: inds(:)
double precision, intent(in) :: src(:)
integer :: bmi_status
end function bmif_set_value_at_indices_double
! Get number of dimensions of the computational grid.
function bmif_get_grid_rank(this, grid, rank) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, intent(out) :: rank
integer :: bmi_status
end function bmif_get_grid_rank
! Get the total number of elements in the computational grid.
function bmif_get_grid_size(this, grid, size) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, intent(out) :: size
integer :: bmi_status
end function bmif_get_grid_size
! Get the grid type as a string.
function bmif_get_grid_type(this, grid, type) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
character(len=*), intent(out) :: type
integer :: bmi_status
end function bmif_get_grid_type
! Get the dimensions of the computational grid.
function bmif_get_grid_shape(this, grid, shape) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, dimension(:), intent(out) :: shape
integer :: bmi_status
end function bmif_get_grid_shape
! Get distance between nodes of the computational grid.
function bmif_get_grid_spacing(this, grid, spacing) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
double precision, dimension(:), intent(out) :: spacing
integer :: bmi_status
end function bmif_get_grid_spacing
! Get coordinates of the origin of the computational grid.
function bmif_get_grid_origin(this, grid, origin) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
double precision, dimension(:), intent(out) :: origin
integer :: bmi_status
end function bmif_get_grid_origin
! Get the x-coordinates of the nodes of a computational grid.
function bmif_get_grid_x(this, grid, x) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
double precision, dimension(:), intent(out) :: x
integer :: bmi_status
end function bmif_get_grid_x
! Get the y-coordinates of the nodes of a computational grid.
function bmif_get_grid_y(this, grid, y) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
double precision, dimension(:), intent(out) :: y
integer :: bmi_status
end function bmif_get_grid_y
! Get the z-coordinates of the nodes of a computational grid.
function bmif_get_grid_z(this, grid, z) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
double precision, dimension(:), intent(out) :: z
integer :: bmi_status
end function bmif_get_grid_z
! Get the number of nodes in an unstructured grid.
function bmif_get_grid_node_count(this, grid, count) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, intent(out) :: count
integer :: bmi_status
end function bmif_get_grid_node_count
! Get the number of edges in an unstructured grid.
function bmif_get_grid_edge_count(this, grid, count) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, intent(out) :: count
integer :: bmi_status
end function bmif_get_grid_edge_count
! Get the number of faces in an unstructured grid.
function bmif_get_grid_face_count(this, grid, count) result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, intent(out) :: count
integer :: bmi_status
end function bmif_get_grid_face_count
! Get the edge-node connectivity.
function bmif_get_grid_edge_nodes(this, grid, edge_nodes) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, dimension(:), intent(out) :: edge_nodes
integer :: bmi_status
end function bmif_get_grid_edge_nodes
! Get the face-edge connectivity.
function bmif_get_grid_face_edges(this, grid, face_edges) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, dimension(:), intent(out) :: face_edges
integer :: bmi_status
end function bmif_get_grid_face_edges
! Get the face-node connectivity.
function bmif_get_grid_face_nodes(this, grid, face_nodes) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, dimension(:), intent(out) :: face_nodes
integer :: bmi_status
end function bmif_get_grid_face_nodes
! Get the number of nodes for each face.
function bmif_get_grid_nodes_per_face(this, grid, nodes_per_face) &
result(bmi_status)
import :: bmi
class(bmi), intent(in) :: this
integer, intent(in) :: grid
integer, dimension(:), intent(out) :: nodes_per_face
integer :: bmi_status
end function bmif_get_grid_nodes_per_face
end interface
end module bmif_2_0