1!
2! CDDL HEADER START
3!
4! The contents of this file are subject to the terms of the Common Development
5! and Distribution License Version 1.0 (the "License").
6!
7! You can obtain a copy of the license at
8! http://www.opensource.org/licenses/CDDL-1.0.  See the License for the
9! specific language governing permissions and limitations under the License.
10!
11! When distributing Covered Code, include this CDDL HEADER in each file and
12! include the License file in a prominent location with the name LICENSE.CDDL.
13! If applicable, add the following below this CDDL HEADER, with the fields
14! enclosed by brackets "[]" replaced with your own identifying information:
15!
16! Portions Copyright (c) [yyyy] [name of copyright owner]. All rights reserved.
17!
18! CDDL HEADER END
19!
20
21!
22! Copyright (c) 2012, Regents of the University of Minnesota.  All rights reserved.
23!
24! Contributors:
25!    Valeriu Smirichinski
26!    Ryan S. Elliott
27!    Ellad B. Tadmor
28!
29
30!
31! Release: This file is part of the openkim-api-v1.1.1 package.
32!
33
34
35#include "KIM_API_status.h"
36
37module kim_kinds
38  implicit none
39#ifdef SYSTEM32
40  integer, parameter :: kim_intptr = 4
41#else
42  integer, parameter :: kim_intptr = 8
43#endif
44end module kim_kinds
45
46module kim_api
47  use kim_kinds
48    implicit none
49    integer,parameter :: KIM_KEY_STRING_LENGTH = 64 !
50
51    interface
52        !explicite inteface to C-side code
53        integer function kim_api_init(kimmdl,testname,mdlname)
54          use kim_kinds
55            integer(kind=kim_intptr) :: kimmdl,testname,mdlname
56        end function kim_api_init
57
58        integer function kim_api_model_info(kimmdl,mdlname)
59          use kim_kinds
60            integer(kind=kim_intptr) :: kimmdl,mdlname
61        end function kim_api_model_info
62
63        integer function kim_api_string_init(kimmdl,testname,mdlname)
64          use kim_kinds
65            integer(kind=kim_intptr) :: kimmdl,testname,mdlname
66        end function kim_api_string_init
67
68        function kim_api_get_status_msg_f(errcode)
69          use kim_kinds
70        integer ::errcode
71        integer(kind=kim_intptr) :: kimmdl,kim_api_get_status_msg_f
72        end function kim_api_get_status_msg_f
73
74        function kim_api_get_model_kim_str(modelname,lenstr,ier)
75          use kim_kinds
76        integer(kind=kim_intptr)::kim_api_get_model_kim_str,modelname
77        integer :: ier,lenstr
78        end function kim_api_get_model_kim_str
79
80        function kim_api_get_model_buffer_f(kimmdl,ier)
81          use kim_kinds
82        integer ::ier
83        integer(kind=kim_intptr) :: kimmdl,kim_api_get_model_buffer_f
84        end function kim_api_get_model_buffer_f
85
86
87        function kim_api_get_test_buffer_f(kimmdl,ier)
88          use kim_kinds
89        integer ::ier
90        integer(kind=kim_intptr) :: kimmdl,kim_api_get_test_buffer_f
91        end function kim_api_get_test_buffer_f
92
93
94        function kim_api_is_half_neighbors_f(kimmdl,ier)
95          use kim_kinds
96        integer ::ier
97        integer(kind=kim_intptr) :: kimmdl,kim_api_is_half_neighbors_f
98        end function kim_api_is_half_neighbors_f
99
100
101        subroutine kim_api_set_model_buffer_f(kimmdl,ob,ier)
102          use kim_kinds
103        integer ::ier
104        integer(kind=kim_intptr) :: kimmdl,ob
105        end subroutine kim_api_set_model_buffer_f
106
107        subroutine kim_api_set_test_buffer_f(kimmdl,ob,ier)
108          use kim_kinds
109        integer ::ier
110        integer(kind=kim_intptr) :: kimmdl,ob
111        end subroutine kim_api_set_test_buffer_f
112
113
114        integer function kim_api_process_dedr_f(pkim,de,r,dx,i,j)
115          use kim_kinds
116                integer(kind=kim_intptr) :: pkim,dx
117                real*8 :: de,r
118                integer ::i,j
119        end function kim_api_process_dedr_f
120
121        integer function kim_api_process_d2edr2_f(pkim,de,r,dx,i,j)
122          use kim_kinds
123                integer(kind=kim_intptr) :: pkim,dx,r,i,j
124                real*8 :: de
125        end function kim_api_process_d2edr2_f
126
127
128    integer function kim_api_report_error(ln,fl,usermsg,ier)
129      use kim_kinds
130      integer :: ln,ier
131      integer(kind=kim_intptr)::fl,usermsg
132    end function kim_api_report_error
133
134
135         integer function kim_api_get_model_index_shift_f(kimmdl)
136           use kim_kinds
137            integer(kind=kim_intptr) :: kimmdl
138        end function kim_api_get_model_index_shift_f
139
140
141
142
143
144        integer function kim_api_set_data(kimmdl,nm, size, dt)
145          use kim_kinds
146            integer(kind=kim_intptr) :: kimmdl,nm,  size, dt
147        end function kim_api_set_data
148
149
150
151        function kim_api_get_data(kimmdl,nm,error)
152          use kim_kinds
153            integer(kind=kim_intptr) :: kimmdl,nm,kim_api_get_data
154            integer::error
155        end function kim_api_get_data
156
157
158
159
160#ifdef FORTRAN2003
161        function kim_api_get_data_cptr(kimmdl,nm)
162          use kim_kinds
163            use iso_c_binding
164           type (c_ptr)::kim_api_get_data_cptr
165            integer(kind=kim_intptr) :: kimmdl,nm
166        end function kim_api_get_data_cptr
167#endif
168
169        function kim_api_get_size(kimmdl,nm,error)
170          use kim_kinds
171            integer(kind=kim_intptr) :: kimmdl,nm,kim_api_get_size
172            integer::error
173        end function kim_api_get_size
174
175        function kim_api_get_neigh_mode_f(kimmdl,error)
176          use kim_kinds
177            integer(kind=kim_intptr) :: kimmdl,kim_api_get_neigh_mode_f
178            integer::error
179        end function kim_api_get_neigh_mode_f
180
181        function kim_api_get_shape(kimmdl,nm, shapea,error)
182          use kim_kinds
183            integer(kind=kim_intptr) :: kimmdl,nm, shapea,kim_api_get_shape
184            integer::error
185        end function kim_api_get_shape
186
187        subroutine kim_api_set_shape(kimmdl,nm,shapea,rank,error)
188          use kim_kinds
189        integer(kind=kim_intptr) :: kimmdl,nm,shapea
190        integer::rank,error
191        end subroutine kim_api_set_shape
192
193        function kim_api_get_model_partcl_typs_f(kimmdl,natypes,error)
194          use kim_kinds
195            integer(kind=kim_intptr) :: kimmdl,kim_api_get_model_partcl_typs_f
196            integer::natypes,error
197        end function kim_api_get_model_partcl_typs_f
198
199        function kim_api_get_test_partcl_typs_f(kimmdl,natypes,error)
200          use kim_kinds
201            integer(kind=kim_intptr) :: kimmdl,kim_api_get_test_partcl_typs_f
202            integer::natypes,error
203        end function kim_api_get_test_partcl_typs_f
204
205        function kim_api_get_params_f(kimmdl,nvpar,error)
206          use kim_kinds
207            integer(kind=kim_intptr) :: kimmdl,kim_api_get_params_f
208            integer::nvpar,error
209        end function kim_api_get_params_f
210
211        function kim_api_get_free_params_f(kimmdl,nvpar,error)
212          use kim_kinds
213            integer(kind=kim_intptr) :: kimmdl,kim_api_get_free_params_f
214            integer::nvpar,error
215        end function kim_api_get_free_params_f
216
217        function kim_api_get_fixed_params_f(kimmdl,nvpar,error)
218          use kim_kinds
219            integer(kind=kim_intptr) :: kimmdl,kim_api_get_fixed_params_f
220            integer::nvpar,error
221        end function kim_api_get_fixed_params_f
222
223        function kim_api_get_nbc_method_f(kimmdl,error)
224          use kim_kinds
225            integer(kind=kim_intptr) :: kimmdl,kim_api_get_nbc_method_f
226            integer::error
227        end function kim_api_get_nbc_method_f
228
229        function kim_api_get_partcl_type_code(kimmdl,nm,error)
230          use kim_kinds
231            integer(kind=kim_intptr) :: kimmdl, nm
232            integer::error,kim_api_get_partcl_type_code
233        end function kim_api_get_partcl_type_code
234
235        subroutine kim_api_set_partcl_type_code(kimmdl,nm,code,error)
236          use kim_kinds
237            integer(kind=kim_intptr) :: kimmdl, nm
238            integer::code, error
239        end subroutine kim_api_set_partcl_type_code
240
241
242
243        integer function kim_api_get_neigh_f(kimmdl,mode,request, atom, numnei, pnei1atom, pRij)
244          use kim_kinds
245        integer :: mode,request,atom,numnei
246        integer(kind=kim_intptr) :: kimmdl,pnei1atom,pRij
247        end function kim_api_get_neigh_f
248
249
250        integer function kim_api_get_compute(kimmdl,nm,error)
251          use kim_kinds
252            integer(kind=kim_intptr) :: kimmdl,nm
253            integer::error
254        end function kim_api_get_compute
255
256        integer function kim_api_get_index(kimmdl,nm,error)
257          use kim_kinds
258            integer(kind=kim_intptr) :: kimmdl,nm
259            integer::error
260        end function kim_api_get_index
261
262        function kim_api_get_data_by_index(kimmdl,I,error)
263          use kim_kinds
264            integer(kind=kim_intptr) :: kimmdl,kim_api_get_data_by_index
265            integer :: I,error
266        end function kim_api_get_data_by_index
267
268        function kim_api_get_size_by_index(kimmdl,I,error)
269          use kim_kinds
270            integer(kind=kim_intptr) :: kimmdl,kim_api_get_size_by_index
271            integer :: I,error
272        end function kim_api_get_size_by_index
273
274        function kim_api_get_shape_by_index(kimmdl,I,shapea,error)
275          use kim_kinds
276            integer(kind=kim_intptr) :: kimmdl,shapea,kim_api_get_shape_by_index
277            integer  :: I,error
278        end function kim_api_get_shape_by_index
279
280        integer function kim_api_get_compute_by_index(kimmdl,I,error)
281          use kim_kinds
282            integer(kind=kim_intptr) :: kimmdl
283            integer :: I,error
284        end function kim_api_get_compute_by_index
285
286
287        subroutine kim_api_allocate(kimmdl,natoms,ntypes,error)
288          use kim_kinds
289            integer(kind=kim_intptr) :: kimmdl
290            integer  :: natoms,ntypes,error
291        end subroutine kim_api_allocate
292
293        subroutine kim_api_free(kimmdl,error)
294          use kim_kinds
295            integer(kind=kim_intptr) :: kimmdl
296        integer::error
297        end subroutine kim_api_free
298
299        subroutine kim_api_print(kimmdl,error)
300          use kim_kinds
301            integer(kind=kim_intptr) :: kimmdl
302            integer ::error
303        end subroutine kim_api_print
304
305        integer function kim_api_model_compute_f(kimmdl)
306          use kim_kinds
307            integer(kind=kim_intptr) :: kimmdl
308        end function kim_api_model_compute_f
309
310        integer function kim_api_model_destroy_f(kimmdl)
311          use kim_kinds
312            integer(kind=kim_intptr) :: kimmdl
313        end function kim_api_model_destroy_f
314
315        integer function kim_api_model_reinit_f(kimmdl)
316          use kim_kinds
317            integer(kind=kim_intptr) :: kimmdl
318        end function kim_api_model_reinit_f
319
320        integer function kim_api_model_init_f(kimmdl)
321          use kim_kinds
322                integer(kind=kim_intptr) :: kimmdl
323        end function kim_api_model_init_f
324
325        real*8 function kim_api_get_scale_conversion(u_from,u_to,error)
326          use kim_kinds
327                integer(kind=kim_intptr) :: u_from,u_to
328                integer:: error
329        end function kim_api_get_scale_conversion
330
331        integer function kim_api_get_unit_handling_f(kimmdl,error)
332          use kim_kinds
333                integer(kind=kim_intptr) :: kimmdl
334                integer ::error
335        end function kim_api_get_unit_handling_f
336
337        function kim_api_get_unit_length_f(kimmdl,error)
338          use kim_kinds
339                integer(kind=kim_intptr) :: kimmdl,kim_api_get_unit_length_f
340                integer ::error
341        end function kim_api_get_unit_length_f
342
343        function kim_api_get_unit_energy_f(kimmdl,error)
344          use kim_kinds
345                integer(kind=kim_intptr) :: kimmdl,kim_api_get_unit_energy_f
346                integer ::error
347        end function kim_api_get_unit_energy_f
348
349        function kim_api_get_unit_charge_f(kimmdl,error)
350          use kim_kinds
351                integer(kind=kim_intptr) :: kimmdl,kim_api_get_unit_charge_f
352                integer ::error
353        end function kim_api_get_unit_charge_f
354
355        function kim_api_get_unit_temperature_f(kimmdl,error)
356          use kim_kinds
357                integer(kind=kim_intptr) :: kimmdl,kim_api_get_unit_temperature_f
358                integer ::error
359        end function kim_api_get_unit_temperature_f
360
361        function kim_api_get_unit_time_f(kimmdl,error)
362          use kim_kinds
363                integer(kind=kim_intptr) :: kimmdl,kim_api_get_unit_time_f
364                integer ::error
365        end function kim_api_get_unit_time_f
366
367        real*8 function kim_api_convert_to_act_unit(kimmdl,length,energy,charge,temperature,time, &
368                 length_exponent, energy_exponent, charge_exponent, temperature_exponent, time_exponent, error)
369          use kim_kinds
370             integer(kind=kim_intptr) ::kimmdl,length,energy,charge,temperature,time
371             real*8:: length_exponent, energy_exponent, charge_exponent, temperature_exponent, time_exponent
372             integer::error
373        end function kim_api_convert_to_act_unit
374
375    !subroutines
376
377
378
379
380
381        subroutine kim_api_set_compute(kimmdl,nm,flag,error)
382          use kim_kinds
383            integer(kind=kim_intptr) :: kimmdl,nm
384        integer::flag,error
385        end subroutine kim_api_set_compute
386
387        subroutine kim_api_set_compute_by_index_f(kimmdl,ind,flag,error)
388          use kim_kinds
389            integer(kind=kim_intptr) :: kimmdl
390        integer::ind,flag,error
391        end subroutine kim_api_set_compute_by_index_f
392
393        integer function kim_api_set_data_by_index(kimmdl,I, size,dt)
394          use kim_kinds
395            integer (kind=kim_intptr):: kimmdl,size,dt
396            integer :: I
397        end function kim_api_set_data_by_index
398
399    end interface
400
401 contains
402    !the following subs converts ctypeArray to fortran type array
403    ! ( to pointer to array with descriptor) without 2003 features
404    ! it is analog of c_f_pointer(...)
405    subroutine KIM_to_F90_real_array_2d(ctypeArray,ArrayWithDescriptor,n,m)
406        implicit none
407        integer  :: n,m
408        real*8,target :: ctypeArray(n,m)
409        real*8,pointer ::ArrayWithDescriptor(:,:)
410        ArrayWithDescriptor=>ctypeArray
411    end subroutine KIM_to_F90_real_array_2d
412
413        subroutine KIM_to_F90_real_array_3d(ctypeArray,ArrayWithDescriptor,n,m,l)
414        implicit none
415        integer  :: n,m,l
416        real*8,target :: ctypeArray(n,m,l)
417        real*8,pointer ::ArrayWithDescriptor(:,:,:)
418        ArrayWithDescriptor=>ctypeArray
419    end subroutine KIM_to_F90_real_array_3d
420
421
422    subroutine KIM_to_F90_int_array_2d(ctypeArray,ArrayWithDescriptor,n,m)
423        implicit none
424        integer  :: n,m !! Check if I can make array big like size of integer*8
425        integer,target :: ctypeArray(n,m)
426        integer,pointer ::ArrayWithDescriptor(:,:)
427        ArrayWithDescriptor=>ctypeArray
428    end subroutine KIM_to_F90_int_array_2d
429
430        subroutine KIM_to_F90_int_array_3d(ctypeArray,ArrayWithDescriptor,n,m,l)
431        implicit none
432        integer  :: n,m,l !! Check if I can make array big like size of integer*8
433        integer,target :: ctypeArray(n,m,l)
434        integer,pointer ::ArrayWithDescriptor(:,:,:)
435        ArrayWithDescriptor=>ctypeArray
436    end subroutine KIM_to_F90_int_array_3d
437
438    subroutine KIM_to_F90_int_array_1d(ctypeArray,ArrayWithDescriptor,n)
439        implicit none
440        integer :: n
441        integer,target :: ctypeArray(n)
442        integer,pointer ::ArrayWithDescriptor(:)
443        ArrayWithDescriptor=>ctypeArray
444    end subroutine KIM_to_F90_int_array_1d
445
446    subroutine KIM_to_F90_real_array_1d(ctypeArray,ArrayWithDescriptor,n)
447        implicit none
448        integer :: n
449        real*8,target :: ctypeArray(n)
450        real*8,pointer ::ArrayWithDescriptor(:)
451        ArrayWithDescriptor=>ctypeArray
452    end subroutine KIM_to_F90_real_array_1d
453
454    character (len=KIM_KEY_STRING_LENGTH) function attachnull(str)
455        character (LEN=*) :: str
456        attachnull = str//CHAR(0)
457    end function attachnull
458    integer(kind=kim_intptr) function attachnull_c(str)
459        character (LEN=*) :: str
460
461  ! print*,"proba02"
462        str = trim(str)//CHAR(0)
463  ! print*,"proba03"
464        attachnull_c = loc(str)
465    end function attachnull_c
466
467    character (len=KIM_KEY_STRING_LENGTH) function cstring2fortran(pointercstring)
468        character (len=KIM_KEY_STRING_LENGTH) cstring ; pointer(pointercstring,cstring)
469        cstring2fortran = cstring
470    end function cstring2fortran
471
472    ! fortran 90 wraper for KIM_API interface using cray pointers
473    ! takes care of character string (null atachment) and dope vector for arrays
474    !****************** f90 wraper ********************
475    integer function kim_api_init_f(kimmdl,testname,mdlname)
476            character (len=*) :: testname,mdlname
477            integer(kind=kim_intptr) :: kimmdl
478            character (len=KIM_KEY_STRING_LENGTH)::testnamesnd,mdlnamesnd
479            character (len=KIM_KEY_STRING_LENGTH):: s1,s2
480            pointer(ps1,s1);pointer(ps2,s2)
481            testnamesnd=attachnull(trim(testname))
482            mdlnamesnd=attachnull(trim(mdlname))
483            ps1=loc(testnamesnd)
484            ps2=loc(mdlnamesnd)
485            kim_api_init_f =kim_api_init(kimmdl,ps1,ps2)
486     end function kim_api_init_f
487
488    integer function kim_api_model_info_f(kimmdl,mdlname)
489            character (len=*) :: mdlname
490            integer(kind=kim_intptr) :: kimmdl
491            character (len=KIM_KEY_STRING_LENGTH)::mdlnamesnd
492            character (len=KIM_KEY_STRING_LENGTH):: s2
493            pointer(ps2,s2)
494            mdlnamesnd=attachnull(trim(mdlname))
495            ps2=loc(mdlnamesnd)
496            kim_api_model_info_f =kim_api_model_info(kimmdl,ps2)
497     end function kim_api_model_info_f
498
499    integer function kim_api_string_init_f(kimmdl,testname,mdlname)
500            character (len=*) :: testname,mdlname
501            integer(kind=kim_intptr) :: kimmdl
502            character (len=KIM_KEY_STRING_LENGTH)::mdlnamesnd
503            character (len=KIM_KEY_STRING_LENGTH):: s1,s2
504            pointer(ps1,s1);pointer(ps2,s2)
505            !testnamesnd=attachnull(trim(testname))
506
507            mdlnamesnd=attachnull(trim(mdlname))
508            ps1=loc(testname)
509            ps2=loc(mdlnamesnd)
510            kim_api_string_init_f =kim_api_string_init(kimmdl,ps1,ps2)
511     end function kim_api_string_init_f
512
513
514
515
516
517        integer function kim_api_set_data_f(kimmdl,nm, size, dt)
518            ! returns 1 (true) if successfull
519            ! dt is address (cray pointer to acctual data
520            integer(kind=kim_intptr) :: kimmdl,  size, dt
521            character (len=*) ::nm
522
523            character (len=KIM_KEY_STRING_LENGTH) :: str2send
524            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
525            str2send = attachnull(trim(nm))
526            pstr = loc(str2send)
527
528            kim_api_set_data_f = kim_api_set_data(kimmdl,pstr,size,dt)
529        end function kim_api_set_data_f
530
531
532        integer(kind=kim_intptr) function kim_api_get_data_f(kimmdl,nm,error)
533            integer(kind=kim_intptr) :: kimmdl
534            character (len=*) ::nm
535            character (len=KIM_KEY_STRING_LENGTH) :: str2send
536            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
537            integer::error
538            str2send = attachnull(trim(nm))
539            pstr = loc(str2send)
540            kim_api_get_data_f = kim_api_get_data(kimmdl,pstr,error)
541        end function kim_api_get_data_f
542
543        real*8 function kim_api_get_scale_conversion_f(from,to,error)
544                integer:: error
545                character (len=*) ::from,to
546                character (len=KIM_KEY_STRING_LENGTH) ::sfrom,sto
547                sfrom = attachnull(trim(from))
548                sto = attachnull(trim(to))
549                kim_api_get_scale_conversion_f = kim_api_get_scale_conversion(loc(sfrom),loc(sto),error)
550        end function kim_api_get_scale_conversion_f
551
552        real*8 function kim_api_convert_to_act_unit_f(kimmdl,length,energy,charge,temperature,time, &
553                 length_exponent, energy_exponent, charge_exponent, temperature_exponent, time_exponent, error)
554                integer::error
555                integer(kind=kim_intptr) :: kimmdl
556                real*8 ::length_exponent, energy_exponent, charge_exponent, temperature_exponent, time_exponent
557                character (len=*) ::length,energy,charge,temperature,time
558                character (len=KIM_KEY_STRING_LENGTH) ::slength,senergy,scharge,stemperature,stime
559                slength = attachnull(trim(length))
560                senergy = attachnull(trim(energy))
561                scharge = attachnull(trim(charge))
562                stemperature = attachnull(trim(temperature))
563                stime = attachnull(trim(time))
564
565                kim_api_convert_to_act_unit_f =kim_api_convert_to_act_unit(kimmdl,loc(slength),loc(senergy),loc(scharge),&
566                loc(stemperature), loc(stime),length_exponent, energy_exponent, charge_exponent, &
567                   temperature_exponent, time_exponent, error)
568        end function kim_api_convert_to_act_unit_f
569
570        integer(kind=kim_intptr) function kim_api_get_model_kim_str_f(nm,lenstr,error)
571            character (len=*) ::nm
572            character (len=KIM_KEY_STRING_LENGTH) :: str2send
573            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
574            integer::error,lenstr
575            str2send = attachnull(trim(nm))
576            pstr = loc(str2send)
577            kim_api_get_model_kim_str_f = kim_api_get_model_kim_str(pstr,lenstr,error)
578        end function kim_api_get_model_kim_str_f
579
580        integer function kim_api_get_partcl_type_code_f(kimmdl,nm,error)
581                integer(kind=kim_intptr) :: kimmdl
582                character (len=*) ::nm
583                character (len=KIM_KEY_STRING_LENGTH) :: str2send
584                character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
585                integer::error
586                str2send = attachnull(trim(nm))
587                pstr = loc(str2send)
588                kim_api_get_partcl_type_code_f = kim_api_get_partcl_type_code(kimmdl,pstr,error)
589        end function kim_api_get_partcl_type_code_f
590
591        subroutine kim_api_set_partcl_type_code_f(kimmdl,nm,code,error)
592                integer(kind=kim_intptr) :: kimmdl
593                character (len=*) ::nm
594                character (len=KIM_KEY_STRING_LENGTH) :: str2send
595                character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
596                integer::code, error
597                str2send = attachnull(trim(nm))
598                pstr = loc(str2send)
599                call kim_api_set_partcl_type_code(kimmdl,pstr,code,error)
600        end subroutine kim_api_set_partcl_type_code_f
601
602
603
604#ifdef FORTRAN2003
605        function kim_api_get_data_fc(kimmdl,nm)
606            use iso_c_binding
607            integer(kind=kim_intptr) :: kimmdl
608            type(c_ptr) :: kim_api_get_data_fc
609            character (len=*) ::nm
610            character (len=KIM_KEY_STRING_LENGTH) :: str2send
611            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
612            str2send = attachnull(trim(nm))
613            pstr = loc(str2send)
614            kim_api_get_data_fc = kim_api_get_data_cptr(kimmdl,pstr)
615        end function kim_api_get_data_fc
616#endif
617
618        integer(kind=kim_intptr) function kim_api_get_size_f(kimmdl,nm,error)
619            integer(kind=kim_intptr) :: kimmdl
620            character (len=*) ::nm
621
622            character (len=KIM_KEY_STRING_LENGTH) :: str2send
623            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
624            integer::error
625            str2send = attachnull(trim(nm))
626            pstr = loc(str2send)
627
628            kim_api_get_size_f = kim_api_get_size(kimmdl,pstr,error)
629        end function kim_api_get_size_f
630
631        integer(kind=kim_intptr) function kim_api_get_shape_f(kimmdl,nm, shape,error)
632            ! returns rank and place correct shape (shape has to be allocated )
633            ! if unssaccesfull returns -1
634            integer(kind=kim_intptr) :: kimmdl
635            character (len=*) :: nm
636
637            character (len=KIM_KEY_STRING_LENGTH) :: str2send
638            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
639
640            integer,pointer :: shape(:)
641            integer :: shp(size(shape))
642            integer :: shpst(1); pointer(pshpst,shpst)
643            integer :: rank,i,error
644            str2send = attachnull(trim(nm))
645            pstr = loc(str2send)
646            pshpst = loc(shp)
647            rank = kim_api_get_shape(kimmdl,pstr,pshpst,error)
648            if(rank .eq. 0) then
649                kim_api_get_shape_f=0
650            else if(rank.eq.1) then
651                kim_api_get_shape_f=1
652                shape(1)=shp(1)
653            else if(rank.gt.1) then
654                kim_api_get_shape_f=rank
655                do i=1, rank
656                    shape(i)=shp(rank - i + 1)
657                end do
658            else
659                kim_api_get_shape_f=-1
660            end if
661        end function kim_api_get_shape_f
662
663        subroutine kim_api_set_shape_f(kimmdl,nm, shapef,rankf,error)
664            ! set rank, shape and size of the data in KIM API object
665            ! error=1 if successful , error <= 0 if unsuccessful
666            integer(kind=kim_intptr) :: kimmdl
667            character (len=*) :: nm
668
669            character (len=KIM_KEY_STRING_LENGTH) :: str2send
670            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
671
672            !integer,pointer :: shape(:)
673            integer ::rankf,shapef(rankf)
674            integer :: shpst(1); pointer(pshpst,shpst)
675            integer :: i,error
676            integer :: shp(rankf)
677            str2send = attachnull(trim(nm))
678            pstr = loc(str2send)
679            pshpst = loc(shp)
680            if(rankf .gt. 1) then
681                do i=1, rankf
682                    shp(i)=shapef(rankf - i + 1) !transpose shape
683                end do
684            end if
685            call kim_api_set_shape(kimmdl,pstr,pshpst,rankf,error)
686        end subroutine kim_api_set_shape_f
687
688        integer function kim_api_get_compute_f(kimmdl,nm,error)
689            integer(kind=kim_intptr) :: kimmdl
690            character (len=*) ::nm
691            integer::error
692            character (len=KIM_KEY_STRING_LENGTH) :: str2send
693            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
694            str2send = attachnull(trim(nm))
695            pstr = loc(str2send)
696
697            kim_api_get_compute_f =kim_api_get_compute(kimmdl,pstr,error)
698        end function kim_api_get_compute_f
699
700        integer function kim_api_get_index_f(kimmdl,nm,error)
701            integer(kind=kim_intptr) :: kimmdl
702            character (len=*) ::nm
703            character (len=KIM_KEY_STRING_LENGTH) :: str2send
704            character (len=KIM_KEY_STRING_LENGTH) :: str; pointer (pstr,str)
705            integer::error
706            str2send = attachnull(trim(nm))
707            pstr = loc(str2send)
708            kim_api_get_index_f = kim_api_get_index(kimmdl,pstr,error)
709        end function kim_api_get_index_f
710
711
712        integer(kind=kim_intptr) function kim_api_get_data_by_index_f(kimmdl,I,error)
713            integer(kind=kim_intptr) :: kimmdl
714            integer :: I,error
715            kim_api_get_data_by_index_f = kim_api_get_data_by_index(kimmdl,I,error)
716        end function kim_api_get_data_by_index_f
717
718        integer(kind=kim_intptr) function kim_api_get_size_by_index_f(kimmdl,I,error)
719            integer(kind=kim_intptr) :: kimmdl
720            integer :: I,error
721            kim_api_get_size_by_index_f = kim_api_get_size_by_index(kimmdl,I,error)
722        end function kim_api_get_size_by_index_f
723
724        integer(kind=kim_intptr) function kim_api_get_shape_by_index_f(kimmdl,I,shape,error)
725            integer(kind=kim_intptr) :: kimmdl
726            integer  :: I,error
727            integer,pointer :: shape(:)
728            integer :: shp(size(shape))
729            integer :: shpst(1); pointer(pshpst,shpst)
730            integer :: rank,ii
731            pshpst = loc(shp)
732            rank = kim_api_get_shape_by_index(kimmdl,I,pshpst,error)
733            if(rank .eq. 0) then
734                kim_api_get_shape_by_index_f=0
735            else if(rank.eq.1) then
736                kim_api_get_shape_by_index_f=1
737                shape(1)=shp(1)
738            else if(rank.gt.1) then
739                kim_api_get_shape_by_index_f=rank
740                do ii=1, rank
741                    shape(ii)=shp(rank - ii + 1)
742                end do
743            else
744                kim_api_get_shape_by_index_f=-1
745            end if
746        end function kim_api_get_shape_by_index_f
747
748        integer function kim_api_get_compute_by_index_f(kimmdl,I,error)
749            integer(kind=kim_intptr) :: kimmdl
750            integer :: I,error
751            kim_api_get_compute_by_index_f = kim_api_get_compute_by_index(kimmdl,I,error)
752        end function kim_api_get_compute_by_index_f
753
754        !subroutines
755        subroutine kim_api_allocate_f(kimmdl,natoms,ntypes,error)
756            integer(kind=kim_intptr) :: kimmdl
757            integer  :: natoms,ntypes,error
758            call kim_api_allocate(kimmdl,natoms,ntypes,error)
759        end subroutine kim_api_allocate_f
760
761        subroutine kim_api_free_f(kimmdl,error)
762            integer(kind=kim_intptr) :: kimmdl
763            integer::error
764            call kim_api_free(kimmdl,error)
765        end subroutine kim_api_free_f
766
767        subroutine kim_api_print_f(kimmdl,error)
768            integer(kind=kim_intptr) :: kimmdl
769            integer::error
770            call kim_api_print(kimmdl,error)
771        end subroutine kim_api_print_f
772    !subroutine kim_api_model_calculate(subroutine *kimmdl) //method is not implemented yet
773
774
775
776
777        integer function kim_api_report_error_f(ln,fl,usermsg,ier)
778            integer:: ln,ier
779            character(len=*)::fl,usermsg
780            character(len=128)::fltmp,umsgtmp
781            fltmp=fl//CHAR(0)
782            umsgtmp=usermsg//char(0)
783            kim_api_report_error_f = kim_api_report_error(ln,loc(fltmp),loc(umsgtmp),ier)
784
785        end function kim_api_report_error_f
786
787        subroutine kim_api_set_compute_f(kimmdl,nm,flag,error)
788            integer(kind=kim_intptr) :: kimmdl
789            character (len=*) :: nm
790            integer :: flag, error
791            character (len=KIM_KEY_STRING_LENGTH) ::us
792            character :: str(1); pointer(pstr,str)
793            us = attachnull(nm)
794            pstr =loc(us)
795
796            call kim_api_set_compute(kimmdl,pstr,flag,error)
797        end subroutine kim_api_set_compute_f
798
799        integer function kim_api_set_data_by_index_f(kimmdl,I, size,dt)
800            integer(kind=kim_intptr) :: kimmdl,size,dt
801            integer :: I
802            kim_api_set_data_by_index_f = kim_api_set_data_by_index(kimmdl,I, size,dt)
803        end function kim_api_set_data_by_index_f
804
805        subroutine kim_api_setm_data_f(kimmdl,error, nm1,sz1,dt1,k1, nm2,sz2,dt2,k2, nm3,sz3,dt3,k3, nm4,sz4,dt4,k4,&
806         nm5,sz5,dt5,k5,   nm6,sz6,dt6,k6,  nm7,sz7,dt7,k7, nm8,sz8,dt8,k8, nm9,sz9,dt9,k9, nm10,sz10,dt10,k10,&
807         nm11,sz11,dt11,k11, nm12,sz12,dt12,k12, nm13,sz13,dt13,k13,  nm14,sz14,dt14,k14, nm15,sz15,dt15,k15)
808            integer(kind=kim_intptr) :: kimmdl;  integer error
809
810            character(len=40) ::msg="kim_api_setm_data_f"
811
812            !    declare optional arguments
813            character(len=*)          :: nm1;   integer(kind=kim_intptr)          :: sz1, dt1
814            character(len=*),optional:: nm2;   integer(kind=kim_intptr),optional:: sz2, dt2
815            character(len=*),optional:: nm3;   integer(kind=kim_intptr),optional:: sz3, dt3
816            character(len=*),optional:: nm4;   integer(kind=kim_intptr),optional:: sz4, dt4
817            character(len=*),optional:: nm5;   integer(kind=kim_intptr),optional:: sz5, dt5
818            character(len=*),optional:: nm6;   integer(kind=kim_intptr),optional:: sz6, dt6
819            character(len=*),optional:: nm7;   integer(kind=kim_intptr),optional:: sz7, dt7
820            character(len=*),optional:: nm8;   integer(kind=kim_intptr),optional:: sz8, dt8
821            character(len=*),optional:: nm9;   integer(kind=kim_intptr),optional:: sz9, dt9
822            character(len=*),optional:: nm10;  integer(kind=kim_intptr),optional:: sz10,dt10
823            character(len=*),optional:: nm11;  integer(kind=kim_intptr),optional:: sz11,dt11
824            character(len=*),optional:: nm12;  integer(kind=kim_intptr),optional:: sz12,dt12
825            character(len=*),optional:: nm13;  integer(kind=kim_intptr),optional:: sz13,dt13
826            character(len=*),optional:: nm14;  integer(kind=kim_intptr),optional:: sz14,dt14
827            character(len=*),optional:: nm15;  integer(kind=kim_intptr),optional:: sz15,dt15
828            integer ::k1;
829            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
830            !
831            if(k1.ne.0.and.k1.ne.1) then
832               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
833               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1, nm1 ).lt.KIM_STATUS_OK) return
834            endif
835            if(k1.eq.1) then
836               error = kim_api_set_data_f(kimmdl,nm1,sz1,dt1);
837               if (errcheck_mltpl(error,msg,1, nm1).lt.KIM_STATUS_OK) return
838            endif
839
840            !check rest of the arguments
841            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
842            if( present(nm2 ).and.(.not.present(sz2 ).or..not.present(dt2 ).or..not.present(k2))) then
843                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
844            elseif( present(nm3 ).and.(.not.present(sz3 ).or..not.present(dt3 ).or..not.present(k3))) then
845                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
846            elseif( present(nm4 ).and.(.not.present(sz4 ).or..not.present(dt4 ).or..not.present(k4))) then
847                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
848            elseif( present(nm5 ).and.(.not.present(sz5 ).or..not.present(dt5 ).or..not.present(k5))) then
849                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
850            elseif( present(nm6 ).and.(.not.present(sz6 ).or..not.present(dt6 ).or..not.present(k6))) then
851                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
852            elseif( present(nm7 ).and.(.not.present(sz7 ).or..not.present(dt7 ).or..not.present(k7))) then
853                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
854            elseif( present(nm8 ).and.(.not.present(sz8 ).or..not.present(dt8 ).or..not.present(k8))) then
855                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
856            elseif( present(nm9 ).and.(.not.present(sz9 ).or..not.present(dt9 ).or..not.present(k9))) then
857                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
858            elseif( present(nm10 ).and.(.not.present(sz10 ).or..not.present(dt10 ).or..not.present(k10))) then
859                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
860            elseif( present(nm11 ).and.(.not.present(sz11 ).or..not.present(dt11 ).or..not.present(k11))) then
861                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
862            elseif( present(nm12 ).and.(.not.present(sz12 ).or..not.present(dt12 ).or..not.present(k12))) then
863                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
864            elseif( present(nm13 ).and.(.not.present(sz13 ).or..not.present(dt13 ).or..not.present(k13))) then
865                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
866            elseif( present(nm14 ).and.(.not.present(sz14 ).or..not.present(dt14 ).or..not.present(k14))) then
867                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
868            elseif( present(nm15 ).and.(.not.present(sz15 ).or..not.present(dt15 ).or..not.present(k15))) then
869                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
870            endif
871
872            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
873            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
874                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
875            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
876                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
877            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
878                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
879            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
880                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
881            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
882                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
883            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
884                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
885            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
886                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
887            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
888                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
889            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
890                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
891            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
892                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
893            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
894                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
895            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
896                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
897            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
898                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
899            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
900                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
901            endif
902
903            !process arguments
904            error = KIM_STATUS_OK
905 if(present(nm2).and.k2.eq.1) error=kim_api_set_data_f(kimmdl,nm2,sz2,dt2);
906            if(errcheck_mltpl(error,msg,2,nm2).lt.KIM_STATUS_OK) return
907 if(present(nm3).and.k3.eq.1) error=kim_api_set_data_f(kimmdl,nm3,sz3,dt3);
908            if(errcheck_mltpl(error,msg,3,nm3).lt.KIM_STATUS_OK) return
909 if(present(nm4).and.k4.eq.1) error=kim_api_set_data_f(kimmdl,nm4,sz4,dt4);
910            if(errcheck_mltpl(error,msg,4,nm4).lt.KIM_STATUS_OK) return
911 if(present(nm5).and.k5.eq.1) error=kim_api_set_data_f(kimmdl,nm5,sz5,dt5);
912            if(errcheck_mltpl(error,msg,5,nm5).lt.KIM_STATUS_OK) return
913 if(present(nm6).and.k6.eq.1) error=kim_api_set_data_f(kimmdl,nm6,sz6,dt6);
914            if(errcheck_mltpl(error,msg,6,nm6).lt.KIM_STATUS_OK) return
915 if(present(nm7).and.k7.eq.1) error=kim_api_set_data_f(kimmdl,nm7,sz7,dt7);
916            if(errcheck_mltpl(error,msg,7,nm7).lt.KIM_STATUS_OK) return
917 if(present(nm8).and.k8.eq.1) error=kim_api_set_data_f(kimmdl,nm8,sz8,dt8);
918            if(errcheck_mltpl(error,msg,8,nm8).lt.KIM_STATUS_OK) return
919 if(present(nm9).and.k9.eq.1) error=kim_api_set_data_f(kimmdl,nm9,sz9,dt9);
920            if(errcheck_mltpl(error,msg,9,nm9).lt.KIM_STATUS_OK) return
921 if(present(nm10).and.k10.eq.1) error=kim_api_set_data_f(kimmdl,nm10,sz10,dt10);
922            if(errcheck_mltpl(error,msg,10,nm10).lt.KIM_STATUS_OK) return
923 if(present(nm11).and.k11.eq.1) error=kim_api_set_data_f(kimmdl,nm11,sz11,dt11);
924            if(errcheck_mltpl(error,msg,11,nm11).lt.KIM_STATUS_OK) return
925 if(present(nm12).and.k12.eq.1) error=kim_api_set_data_f(kimmdl,nm12,sz12,dt12);
926            if(errcheck_mltpl(error,msg,12,nm12).lt.KIM_STATUS_OK) return
927 if(present(nm13).and.k13.eq.1) error=kim_api_set_data_f(kimmdl,nm13,sz13,dt13);
928            if(errcheck_mltpl(error,msg,13,nm13).lt.KIM_STATUS_OK) return
929 if(present(nm14).and.k14.eq.1) error=kim_api_set_data_f(kimmdl,nm14,sz14,dt14);
930            if(errcheck_mltpl(error,msg,14,nm14).lt.KIM_STATUS_OK) return
931 if(present(nm15).and.k15.eq.1) error=kim_api_set_data_f(kimmdl,nm15,sz15,dt15);
932            if(errcheck_mltpl(error,msg,15,nm15).lt.KIM_STATUS_OK) return
933        end subroutine kim_api_setm_data_f
934
935        subroutine kim_api_setm_data_by_index_f(kimmdl,error,nm1,sz1,dt1,k1, nm2,sz2,dt2,k2, nm3,sz3,dt3,k3, nm4,sz4,dt4,k4,&
936         nm5,sz5,dt5,k5,   nm6,sz6,dt6,k6,  nm7,sz7,dt7,k7, nm8,sz8,dt8,k8, nm9,sz9,dt9,k9, nm10,sz10,dt10,k10,&
937         nm11,sz11,dt11,k11, nm12,sz12,dt12,k12, nm13,sz13,dt13,k13,  nm14,sz14,dt14,k14, nm15,sz15,dt15,k15)
938
939            integer(kind=kim_intptr) :: kimmdl;  integer error
940
941            character(len=40) ::msg="kim_api_setm_data_by_index_f"
942
943            !    declare optional arguments
944            integer          :: nm1;   integer(kind=kim_intptr)          :: sz1, dt1
945            integer,optional:: nm2;   integer(kind=kim_intptr),optional:: sz2, dt2
946            integer,optional:: nm3;   integer(kind=kim_intptr),optional:: sz3, dt3
947            integer,optional:: nm4;   integer(kind=kim_intptr),optional:: sz4, dt4
948            integer,optional:: nm5;   integer(kind=kim_intptr),optional:: sz5, dt5
949            integer,optional:: nm6;   integer(kind=kim_intptr),optional:: sz6, dt6
950            integer,optional:: nm7;   integer(kind=kim_intptr),optional:: sz7, dt7
951            integer,optional:: nm8;   integer(kind=kim_intptr),optional:: sz8, dt8
952            integer,optional:: nm9;   integer(kind=kim_intptr),optional:: sz9, dt9
953            integer,optional:: nm10;  integer(kind=kim_intptr),optional:: sz10,dt10
954            integer,optional:: nm11;  integer(kind=kim_intptr),optional:: sz11,dt11
955            integer,optional:: nm12;  integer(kind=kim_intptr),optional:: sz12,dt12
956            integer,optional:: nm13;  integer(kind=kim_intptr),optional:: sz13,dt13
957            integer,optional:: nm14;  integer(kind=kim_intptr),optional:: sz14,dt14
958            integer,optional:: nm15;  integer(kind=kim_intptr),optional:: sz15,dt15
959            integer ::k1;
960            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
961            !
962            if(k1.ne.0.and.k1.ne.1) then
963               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
964               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1,"",nm1 ).lt.KIM_STATUS_OK) return
965            endif
966            if(k1.eq.1) then
967               error = kim_api_set_data_by_index_f(kimmdl,nm1,sz1,dt1);
968               if (errcheck_mltpl(error,msg,1,"", nm1).lt.KIM_STATUS_OK) return
969            endif
970
971            !check rest of the arguments
972            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
973            if( present(nm2 ).and.(.not.present(sz2 ).or..not.present(dt2 ))) then
974                if(errcheck_mltpl(error,msg,2,"", nm2 ).lt.KIM_STATUS_OK) return
975            elseif( present(nm3 ).and.(.not.present(sz3 ).or..not.present(dt3 ))) then
976                if(errcheck_mltpl(error,msg,3,"", nm3 ).lt.KIM_STATUS_OK) return
977            elseif( present(nm4 ).and.(.not.present(sz4 ).or..not.present(dt4 ))) then
978                if(errcheck_mltpl(error,msg,4,"", nm4 ).lt.KIM_STATUS_OK) return
979            elseif( present(nm5 ).and.(.not.present(sz5 ).or..not.present(dt5 ))) then
980                if(errcheck_mltpl(error,msg,5,"", nm5 ).lt.KIM_STATUS_OK) return
981            elseif( present(nm6 ).and.(.not.present(sz6 ).or..not.present(dt6 ))) then
982                if(errcheck_mltpl(error,msg,6,"", nm6 ).lt.KIM_STATUS_OK) return
983            elseif( present(nm7 ).and.(.not.present(sz7 ).or..not.present(dt7 ))) then
984                if(errcheck_mltpl(error,msg,7,"", nm7 ).lt.KIM_STATUS_OK) return
985            elseif( present(nm8 ).and.(.not.present(sz8 ).or..not.present(dt8 ))) then
986                if(errcheck_mltpl(error,msg,8,"", nm8 ).lt.KIM_STATUS_OK) return
987            elseif( present(nm9 ).and.(.not.present(sz9 ).or..not.present(dt9 ))) then
988                if(errcheck_mltpl(error,msg,9,"", nm9 ).lt.KIM_STATUS_OK) return
989            elseif( present(nm10 ).and.(.not.present(sz10 ).or..not.present(dt10 ))) then
990                if(errcheck_mltpl(error,msg,10,"", nm10 ).lt.KIM_STATUS_OK) return
991            elseif( present(nm11 ).and.(.not.present(sz11 ).or..not.present(dt11 ))) then
992                if(errcheck_mltpl(error,msg,11,"", nm11 ).lt.KIM_STATUS_OK) return
993            elseif( present(nm12 ).and.(.not.present(sz12 ).or..not.present(dt12 ))) then
994                if(errcheck_mltpl(error,msg,12,"", nm12 ).lt.KIM_STATUS_OK) return
995            elseif( present(nm13 ).and.(.not.present(sz13 ).or..not.present(dt13 ))) then
996                if(errcheck_mltpl(error,msg,13,"",nm13 ).lt.KIM_STATUS_OK) return
997            elseif( present(nm14 ).and.(.not.present(sz14 ).or..not.present(dt14 ))) then
998                if(errcheck_mltpl(error,msg,14,"", nm14 ).lt.KIM_STATUS_OK) return
999            elseif( present(nm15 ).and.(.not.present(sz15 ).or..not.present(dt15 ))) then
1000                if(errcheck_mltpl(error,msg,15,"", nm15 ).lt.KIM_STATUS_OK) return
1001            endif
1002
1003            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1004            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1005                if(errcheck_mltpl(error,msg,2,"", nm2 ).lt.KIM_STATUS_OK) return
1006            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1007                if(errcheck_mltpl(error,msg,3,"", nm3 ).lt.KIM_STATUS_OK) return
1008            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1009                if(errcheck_mltpl(error,msg,4,"", nm4 ).lt.KIM_STATUS_OK) return
1010            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1011                if(errcheck_mltpl(error,msg,5,"", nm5 ).lt.KIM_STATUS_OK) return
1012            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1013                if(errcheck_mltpl(error,msg,6,"", nm6 ).lt.KIM_STATUS_OK) return
1014            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1015                if(errcheck_mltpl(error,msg,7,"", nm7 ).lt.KIM_STATUS_OK) return
1016            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1017                if(errcheck_mltpl(error,msg,8,"", nm8 ).lt.KIM_STATUS_OK) return
1018            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1019                if(errcheck_mltpl(error,msg,9,"", nm9 ).lt.KIM_STATUS_OK) return
1020            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1021                if(errcheck_mltpl(error,msg,10,"", nm10 ).lt.KIM_STATUS_OK) return
1022            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1023                if(errcheck_mltpl(error,msg,11,"", nm11 ).lt.KIM_STATUS_OK) return
1024            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1025                if(errcheck_mltpl(error,msg,12,"", nm12 ).lt.KIM_STATUS_OK) return
1026            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1027                if(errcheck_mltpl(error,msg,13,"", nm13 ).lt.KIM_STATUS_OK) return
1028            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1029                if(errcheck_mltpl(error,msg,14,"", nm14 ).lt.KIM_STATUS_OK) return
1030            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1031                if(errcheck_mltpl(error,msg,15,"", nm15 ).lt.KIM_STATUS_OK) return
1032            endif
1033
1034            !process arguments
1035            error=KIM_STATUS_OK
1036 if(present(nm2).and.k2.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm2,sz2,dt2);
1037            if(errcheck_mltpl(error,msg,2,"",nm2).lt.KIM_STATUS_OK) return
1038 if(present(nm3).and.k3.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm3,sz3,dt3);
1039            if(errcheck_mltpl(error,msg,3,"",nm3).lt.KIM_STATUS_OK) return
1040 if(present(nm4).and.k4.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm4,sz4,dt4);
1041            if(errcheck_mltpl(error,msg,4,"",nm4).lt.KIM_STATUS_OK) return
1042 if(present(nm5).and.k5.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm5,sz5,dt5);
1043            if(errcheck_mltpl(error,msg,5,"",nm5).lt.KIM_STATUS_OK) return
1044 if(present(nm6).and.k6.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm6,sz6,dt6);
1045            if(errcheck_mltpl(error,msg,6,"",nm6).lt.KIM_STATUS_OK) return
1046 if(present(nm7).and.k7.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm7,sz7,dt7);
1047            if(errcheck_mltpl(error,msg,7,"",nm7).lt.KIM_STATUS_OK) return
1048 if(present(nm8).and.k8.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm8,sz8,dt8);
1049            if(errcheck_mltpl(error,msg,8,"",nm8).lt.KIM_STATUS_OK) return
1050 if(present(nm9).and.k9.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm9,sz9,dt9);
1051            if(errcheck_mltpl(error,msg,9,"",nm9).lt.KIM_STATUS_OK) return
1052 if(present(nm10).and.k10.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm10,sz10,dt10);
1053            if(errcheck_mltpl(error,msg,10,"",nm10).lt.KIM_STATUS_OK)return
1054 if(present(nm11).and.k11.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm11,sz11,dt11);
1055            if(errcheck_mltpl(error,msg,11,"",nm11).lt.KIM_STATUS_OK)return
1056 if(present(nm12).and.k12.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm12,sz12,dt12);
1057            if(errcheck_mltpl(error,msg,12,"",nm12).lt.KIM_STATUS_OK)return
1058 if(present(nm13).and.k13.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm13,sz13,dt13);
1059            if(errcheck_mltpl(error,msg,13,"",nm13).lt.KIM_STATUS_OK)return
1060 if(present(nm14).and.k14.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm14,sz14,dt14);
1061            if(errcheck_mltpl(error,msg,14,"",nm14).lt.KIM_STATUS_OK)return
1062 if(present(nm15).and.k15.eq.1) error = kim_api_set_data_by_index_f(kimmdl,nm15,sz15,dt15);
1063            if(errcheck_mltpl(error,msg,15,"",nm15).lt.KIM_STATUS_OK)return
1064        end subroutine kim_api_setm_data_by_index_f
1065
1066        subroutine kim_api_getm_data_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1067         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1068         nm14,dt14,k14, nm15,dt15,k15)
1069            integer(kind=kim_intptr) :: kimmdl;  integer error
1070
1071            character(len=40) ::msg="kim_api_getm_data_f"
1072
1073            !    declare optional arguments
1074            character(len=*)          :: nm1;   integer(kind=kim_intptr)          :: dt1
1075            character(len=*),optional:: nm2;   integer(kind=kim_intptr),optional:: dt2
1076            character(len=*),optional:: nm3;   integer(kind=kim_intptr),optional:: dt3
1077            character(len=*),optional:: nm4;   integer(kind=kim_intptr),optional:: dt4
1078            character(len=*),optional:: nm5;   integer(kind=kim_intptr),optional:: dt5
1079            character(len=*),optional:: nm6;   integer(kind=kim_intptr),optional:: dt6
1080            character(len=*),optional:: nm7;   integer(kind=kim_intptr),optional:: dt7
1081            character(len=*),optional:: nm8;   integer(kind=kim_intptr),optional:: dt8
1082            character(len=*),optional:: nm9;   integer(kind=kim_intptr),optional:: dt9
1083            character(len=*),optional:: nm10;  integer(kind=kim_intptr),optional:: dt10
1084            character(len=*),optional:: nm11;  integer(kind=kim_intptr),optional:: dt11
1085            character(len=*),optional:: nm12;  integer(kind=kim_intptr),optional:: dt12
1086            character(len=*),optional:: nm13;  integer(kind=kim_intptr),optional:: dt13
1087            character(len=*),optional:: nm14;  integer(kind=kim_intptr),optional:: dt14
1088            character(len=*),optional:: nm15;  integer(kind=kim_intptr),optional:: dt15
1089            integer ::k1;
1090            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1091            !
1092            if(k1.ne.0.and.k1.ne.1) then
1093               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1094               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1, nm1 ).lt.KIM_STATUS_OK) return
1095            endif
1096            if(k1.eq.1) then
1097                dt1 = kim_api_get_data_f(kimmdl,nm1,error);
1098                if (errcheck_mltpl(error,msg,1, nm1).lt.KIM_STATUS_OK) return
1099            endif
1100
1101            !check rest of the arguments
1102            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1103            if( present(nm2 ).and.(.not.present(dt2 ))) then
1104                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1105            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1106                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1107            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1108                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1109            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1110                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1111            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1112                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1113            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1114                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1115            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1116                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1117            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1118                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1119            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1120                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1121            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1122                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1123            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1124                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1125            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1126                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1127            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1128                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1129            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1130                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1131            endif
1132
1133            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1134            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1135                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1136            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1137                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1138            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1139                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1140            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1141                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1142            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1143                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1144            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1145                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1146            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1147                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1148            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1149                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1150            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1151                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1152            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1153                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1154            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1155                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1156            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1157                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1158            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1159                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1160            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1161                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1162            endif
1163
1164            !process arguments
1165            error=KIM_STATUS_OK
1166 if(present(nm2).and.k2.eq.1) dt2= kim_api_get_data_f(kimmdl,nm2,error);
1167            if(errcheck_mltpl(error,msg,2,nm2).lt.KIM_STATUS_OK) return
1168 if(present(nm3).and.k3.eq.1) dt3= kim_api_get_data_f(kimmdl,nm3,error);
1169            if(errcheck_mltpl(error,msg,3,nm3).lt.KIM_STATUS_OK) return
1170 if(present(nm4).and.k4.eq.1) dt4= kim_api_get_data_f(kimmdl,nm4,error);
1171            if(errcheck_mltpl(error,msg,4,nm4).lt.KIM_STATUS_OK) return
1172 if(present(nm5).and.k5.eq.1) dt5= kim_api_get_data_f(kimmdl,nm5,error);
1173            if(errcheck_mltpl(error,msg,5,nm5).lt.KIM_STATUS_OK) return
1174 if(present(nm6).and.k6.eq.1) dt6= kim_api_get_data_f(kimmdl,nm6,error);
1175            if(errcheck_mltpl(error,msg,6,nm6).lt.KIM_STATUS_OK) return
1176 if(present(nm7).and.k7.eq.1) dt7= kim_api_get_data_f(kimmdl,nm7,error);
1177            if(errcheck_mltpl(error,msg,7,nm7).lt.KIM_STATUS_OK) return
1178 if(present(nm8).and.k8.eq.1) dt8= kim_api_get_data_f(kimmdl,nm8,error);
1179            if(errcheck_mltpl(error,msg,8,nm8).lt.KIM_STATUS_OK) return
1180 if(present(nm9).and.k9.eq.1) dt9= kim_api_get_data_f(kimmdl,nm9,error);
1181            if(errcheck_mltpl(error,msg,9,nm9).lt.KIM_STATUS_OK) return
1182 if(present(nm10).and.k10.eq.1) dt10= kim_api_get_data_f(kimmdl,nm10,error);
1183            if(errcheck_mltpl(error,msg,10,nm10).lt.KIM_STATUS_OK) return
1184 if(present(nm11).and.k11.eq.1) dt11= kim_api_get_data_f(kimmdl,nm11,error);
1185            if(errcheck_mltpl(error,msg,11,nm11).lt.KIM_STATUS_OK) return
1186 if(present(nm12).and.k12.eq.1) dt12= kim_api_get_data_f(kimmdl,nm12,error);
1187            if(errcheck_mltpl(error,msg,12,nm12).lt.KIM_STATUS_OK) return
1188 if(present(nm13).and.k13.eq.1) dt13= kim_api_get_data_f(kimmdl,nm13,error);
1189            if(errcheck_mltpl(error,msg,13,nm13).lt.KIM_STATUS_OK) return
1190 if(present(nm14).and.k14.eq.1) dt14= kim_api_get_data_f(kimmdl,nm14,error);
1191            if(errcheck_mltpl(error,msg,14,nm14).lt.KIM_STATUS_OK) return
1192 if(present(nm15).and.k15.eq.1) dt15= kim_api_get_data_f(kimmdl,nm15,error);
1193            if(errcheck_mltpl(error,msg,15,nm15).lt.KIM_STATUS_OK) return
1194        end subroutine kim_api_getm_data_f
1195
1196        subroutine kim_api_getm_data_by_index_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1197         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1198         nm14,dt14,k14, nm15,dt15,k15)
1199            integer(kind=kim_intptr) :: kimmdl;  integer error
1200
1201            character(len=40) ::msg="kim_api_getm_data_by_index_f"
1202
1203            !    declare optional arguments
1204            integer          :: nm1;   integer(kind=kim_intptr)          :: dt1
1205            integer,optional:: nm2;   integer(kind=kim_intptr),optional:: dt2
1206            integer,optional:: nm3;   integer(kind=kim_intptr),optional:: dt3
1207            integer,optional:: nm4;   integer(kind=kim_intptr),optional:: dt4
1208            integer,optional:: nm5;   integer(kind=kim_intptr),optional:: dt5
1209            integer,optional:: nm6;   integer(kind=kim_intptr),optional:: dt6
1210            integer,optional:: nm7;   integer(kind=kim_intptr),optional:: dt7
1211            integer,optional:: nm8;   integer(kind=kim_intptr),optional:: dt8
1212            integer,optional:: nm9;   integer(kind=kim_intptr),optional:: dt9
1213            integer,optional:: nm10;  integer(kind=kim_intptr),optional:: dt10
1214            integer,optional:: nm11;  integer(kind=kim_intptr),optional:: dt11
1215            integer,optional:: nm12;  integer(kind=kim_intptr),optional:: dt12
1216            integer,optional:: nm13;  integer(kind=kim_intptr),optional:: dt13
1217            integer,optional:: nm14;  integer(kind=kim_intptr),optional:: dt14
1218            integer,optional:: nm15;  integer(kind=kim_intptr),optional:: dt15
1219            integer ::k1;
1220            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1221            !
1222            if(k1.ne.0.and.k1.ne.1) then
1223               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1224              if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1,"", nm1 ).lt.KIM_STATUS_OK) return
1225            endif
1226            if(k1.eq.1) then
1227                dt1 = kim_api_get_data_by_index_f(kimmdl,nm1,error);
1228                if (errcheck_mltpl(error,msg,1, "", nm1).lt.KIM_STATUS_OK) return
1229            endif
1230            !check rest of the arguments
1231            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1232            if( present(nm2 ).and.(.not.present(dt2 ))) then
1233                if(errcheck_mltpl(error,msg,2, "",nm2 ).lt.KIM_STATUS_OK) return
1234            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1235                if(errcheck_mltpl(error,msg,3, "",nm3 ).lt.KIM_STATUS_OK) return
1236            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1237                if(errcheck_mltpl(error,msg,4, "",nm4 ).lt.KIM_STATUS_OK) return
1238            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1239                if(errcheck_mltpl(error,msg,5, "",nm5 ).lt.KIM_STATUS_OK) return
1240            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1241                if(errcheck_mltpl(error,msg,6, "",nm6 ).lt.KIM_STATUS_OK) return
1242            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1243                if(errcheck_mltpl(error,msg,7, "",nm7 ).lt.KIM_STATUS_OK) return
1244            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1245                if(errcheck_mltpl(error,msg,8, "",nm8 ).lt.KIM_STATUS_OK) return
1246            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1247                if(errcheck_mltpl(error,msg,9, "",nm9 ).lt.KIM_STATUS_OK) return
1248            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1249                if(errcheck_mltpl(error,msg,10, "",nm10 ).lt.KIM_STATUS_OK) return
1250            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1251                if(errcheck_mltpl(error,msg,11, "",nm11 ).lt.KIM_STATUS_OK) return
1252            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1253                if(errcheck_mltpl(error,msg,12, "",nm12 ).lt.KIM_STATUS_OK) return
1254            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1255                if(errcheck_mltpl(error,msg,13, "",nm13 ).lt.KIM_STATUS_OK) return
1256            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1257                if(errcheck_mltpl(error,msg,14, "",nm14 ).lt.KIM_STATUS_OK) return
1258            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1259                if(errcheck_mltpl(error,msg,15, "",nm15 ).lt.KIM_STATUS_OK) return
1260            endif
1261
1262            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1263            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1264                if(errcheck_mltpl(error,msg,2,"", nm2 ).lt.KIM_STATUS_OK) return
1265            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1266                if(errcheck_mltpl(error,msg,3,"", nm3 ).lt.KIM_STATUS_OK) return
1267            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1268                if(errcheck_mltpl(error,msg,4,"", nm4 ).lt.KIM_STATUS_OK) return
1269            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1270                if(errcheck_mltpl(error,msg,5,"", nm5 ).lt.KIM_STATUS_OK) return
1271            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1272                if(errcheck_mltpl(error,msg,6,"", nm6 ).lt.KIM_STATUS_OK) return
1273            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1274                if(errcheck_mltpl(error,msg,7,"", nm7 ).lt.KIM_STATUS_OK) return
1275            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1276                if(errcheck_mltpl(error,msg,8,"", nm8 ).lt.KIM_STATUS_OK) return
1277            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1278                if(errcheck_mltpl(error,msg,9,"", nm9 ).lt.KIM_STATUS_OK) return
1279            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1280                if(errcheck_mltpl(error,msg,10,"", nm10 ).lt.KIM_STATUS_OK) return
1281            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1282                if(errcheck_mltpl(error,msg,11,"", nm11 ).lt.KIM_STATUS_OK) return
1283            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1284                if(errcheck_mltpl(error,msg,12,"", nm12 ).lt.KIM_STATUS_OK) return
1285            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1286                if(errcheck_mltpl(error,msg,13,"", nm13 ).lt.KIM_STATUS_OK) return
1287            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1288                if(errcheck_mltpl(error,msg,14,"", nm14 ).lt.KIM_STATUS_OK) return
1289            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1290                if(errcheck_mltpl(error,msg,15,"", nm15 ).lt.KIM_STATUS_OK) return
1291            endif
1292
1293            !process arguments
1294            error =KIM_STATUS_OK
1295 if(present(nm2).and.k2.eq.1) dt2= kim_api_get_data_by_index_f(kimmdl,nm2,error);
1296            if(errcheck_mltpl(error,msg,2,"",nm2).lt.KIM_STATUS_OK) return
1297 if(present(nm3).and.k3.eq.1) dt3= kim_api_get_data_by_index_f(kimmdl,nm3,error);
1298            if(errcheck_mltpl(error,msg,3,"",nm3).lt.KIM_STATUS_OK) return
1299 if(present(nm4).and.k4.eq.1) dt4= kim_api_get_data_by_index_f(kimmdl,nm4,error);
1300            if(errcheck_mltpl(error,msg,4,"",nm4).lt.KIM_STATUS_OK) return
1301 if(present(nm5).and.k5.eq.1) dt5= kim_api_get_data_by_index_f(kimmdl,nm5,error);
1302            if(errcheck_mltpl(error,msg,5,"",nm5).lt.KIM_STATUS_OK) return
1303 if(present(nm6).and.k6.eq.1) dt6= kim_api_get_data_by_index_f(kimmdl,nm6,error);
1304            if(errcheck_mltpl(error,msg,6,"",nm6).lt.KIM_STATUS_OK) return
1305 if(present(nm7).and.k7.eq.1) dt7= kim_api_get_data_by_index_f(kimmdl,nm7,error);
1306            if(errcheck_mltpl(error,msg,7,"",nm7).lt.KIM_STATUS_OK) return
1307 if(present(nm8).and.k8.eq.1) dt8= kim_api_get_data_by_index_f(kimmdl,nm8,error);
1308            if(errcheck_mltpl(error,msg,8,"",nm8).lt.KIM_STATUS_OK) return
1309 if(present(nm9).and.k9.eq.1) dt9= kim_api_get_data_by_index_f(kimmdl,nm9,error);
1310            if(errcheck_mltpl(error,msg,9,"",nm9).lt.KIM_STATUS_OK) return
1311 if(present(nm10).and.k10.eq.1) dt10= kim_api_get_data_by_index_f(kimmdl,nm10,error);
1312            if(errcheck_mltpl(error,msg,10,"",nm10).lt.KIM_STATUS_OK) return
1313 if(present(nm11).and.k11.eq.1) dt11= kim_api_get_data_by_index_f(kimmdl,nm11,error);
1314            if(errcheck_mltpl(error,msg,11,"",nm11).lt.KIM_STATUS_OK) return
1315 if(present(nm12).and.k12.eq.1) dt12= kim_api_get_data_by_index_f(kimmdl,nm12,error);
1316            if(errcheck_mltpl(error,msg,12,"",nm12).lt.KIM_STATUS_OK) return
1317 if(present(nm13).and.k13.eq.1) dt13= kim_api_get_data_by_index_f(kimmdl,nm13,error);
1318            if(errcheck_mltpl(error,msg,13,"",nm13).lt.KIM_STATUS_OK) return
1319 if(present(nm14).and.k14.eq.1) dt14= kim_api_get_data_by_index_f(kimmdl,nm14,error);
1320            if(errcheck_mltpl(error,msg,14,"",nm14).lt.KIM_STATUS_OK) return
1321 if(present(nm15).and.k15.eq.1) dt15= kim_api_get_data_by_index_f(kimmdl,nm15,error);
1322            if(errcheck_mltpl(error,msg,15,"",nm15).lt.KIM_STATUS_OK) return
1323        end subroutine kim_api_getm_data_by_index_f
1324
1325        subroutine kim_api_getm_compute_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1326         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1327         nm14,dt14,k14, nm15,dt15,k15)
1328            integer(kind=kim_intptr) :: kimmdl;  integer error
1329
1330            character(len=40) ::msg="kim_api_getm_compute_f"
1331
1332            !    declare optional arguments
1333            character(len=*)          :: nm1;   integer          :: dt1
1334            character(len=*),optional:: nm2;   integer,optional:: dt2
1335            character(len=*),optional:: nm3;   integer,optional:: dt3
1336            character(len=*),optional:: nm4;   integer,optional:: dt4
1337            character(len=*),optional:: nm5;   integer,optional:: dt5
1338            character(len=*),optional:: nm6;   integer,optional:: dt6
1339            character(len=*),optional:: nm7;   integer,optional:: dt7
1340            character(len=*),optional:: nm8;   integer,optional:: dt8
1341            character(len=*),optional:: nm9;   integer,optional:: dt9
1342            character(len=*),optional:: nm10;  integer,optional:: dt10
1343            character(len=*),optional:: nm11;  integer,optional:: dt11
1344            character(len=*),optional:: nm12;  integer,optional:: dt12
1345            character(len=*),optional:: nm13;  integer,optional:: dt13
1346            character(len=*),optional:: nm14;  integer,optional:: dt14
1347            character(len=*),optional:: nm15;  integer,optional:: dt15
1348            integer ::k1;
1349            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1350           !
1351            if(k1.ne.0.and.k1.ne.1) then
1352               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1353              if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1, nm1 ).lt.KIM_STATUS_OK) return
1354            endif
1355            if(k1.eq.1) then
1356                dt1 = kim_api_get_compute_f(kimmdl,nm1,error);
1357                if (errcheck_mltpl(error,msg,1, nm1).lt.KIM_STATUS_OK) return
1358            endif
1359
1360            !check rest of the arguments
1361            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1362            if( present(nm2 ).and.(.not.present(dt2 ))) then
1363                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1364            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1365                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1366            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1367                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1368            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1369                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1370            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1371                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1372            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1373                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1374            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1375                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1376            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1377                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1378            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1379                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1380            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1381                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1382            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1383                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1384            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1385                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1386            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1387                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1388            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1389                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1390            endif
1391
1392            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1393            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1394                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1395            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1396                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1397            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1398                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1399            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1400                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1401            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1402                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1403            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1404                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1405            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1406                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1407            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1408                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1409            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1410                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1411            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1412                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1413            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1414                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1415            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1416                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1417            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1418                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1419            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1420                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1421            endif
1422
1423            !process arguments
1424            error=KIM_STATUS_OK
1425 if(present(nm2).and.k2.eq.1) dt2= kim_api_get_compute_f(kimmdl,nm2,error);
1426            if(errcheck_mltpl(error,msg,2,nm2).lt.KIM_STATUS_OK) return
1427 if(present(nm3).and.k3.eq.1) dt3= kim_api_get_compute_f(kimmdl,nm3,error);
1428            if(errcheck_mltpl(error,msg,3,nm3).lt.KIM_STATUS_OK) return
1429 if(present(nm4).and.k4.eq.1) dt4= kim_api_get_compute_f(kimmdl,nm4,error);
1430            if(errcheck_mltpl(error,msg,4,nm4).lt.KIM_STATUS_OK) return
1431 if(present(nm5).and.k5.eq.1) dt5= kim_api_get_compute_f(kimmdl,nm5,error);
1432            if(errcheck_mltpl(error,msg,5,nm5).lt.KIM_STATUS_OK) return
1433 if(present(nm6).and.k6.eq.1) dt6= kim_api_get_compute_f(kimmdl,nm6,error);
1434            if(errcheck_mltpl(error,msg,6,nm6).lt.KIM_STATUS_OK) return
1435 if(present(nm7).and.k7.eq.1) dt7= kim_api_get_compute_f(kimmdl,nm7,error);
1436            if(errcheck_mltpl(error,msg,7,nm7).lt.KIM_STATUS_OK) return
1437 if(present(nm8).and.k8.eq.1) dt8= kim_api_get_compute_f(kimmdl,nm8,error);
1438            if(errcheck_mltpl(error,msg,8,nm8).lt.KIM_STATUS_OK) return
1439 if(present(nm9).and.k9.eq.1) dt9= kim_api_get_compute_f(kimmdl,nm9,error);
1440            if(errcheck_mltpl(error,msg,9,nm9).lt.KIM_STATUS_OK) return
1441 if(present(nm10).and.k10.eq.1) dt10= kim_api_get_compute_f(kimmdl,nm10,error);
1442            if(errcheck_mltpl(error,msg,10,nm10).lt.KIM_STATUS_OK) return
1443 if(present(nm11).and.k11.eq.1) dt11= kim_api_get_compute_f(kimmdl,nm11,error);
1444            if(errcheck_mltpl(error,msg,11,nm11).lt.KIM_STATUS_OK) return
1445 if(present(nm12).and.k12.eq.1) dt12= kim_api_get_compute_f(kimmdl,nm12,error);
1446            if(errcheck_mltpl(error,msg,12,nm12).lt.KIM_STATUS_OK) return
1447 if(present(nm13).and.k13.eq.1) dt13= kim_api_get_compute_f(kimmdl,nm13,error);
1448            if(errcheck_mltpl(error,msg,13,nm13).lt.KIM_STATUS_OK) return
1449 if(present(nm14).and.k14.eq.1) dt14= kim_api_get_compute_f(kimmdl,nm14,error);
1450            if(errcheck_mltpl(error,msg,14,nm14).lt.KIM_STATUS_OK) return
1451 if(present(nm15).and.k15.eq.1) dt15= kim_api_get_compute_f(kimmdl,nm15,error);
1452            if(errcheck_mltpl(error,msg,15,nm15).lt.KIM_STATUS_OK) return
1453        end subroutine kim_api_getm_compute_f
1454
1455        subroutine kim_api_getm_compute_by_index_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1456         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1457         nm14,dt14,k14, nm15,dt15,k15)
1458            integer(kind=kim_intptr) :: kimmdl;  integer error
1459
1460            character(len=40) ::msg="kim_api_getm_compute_by_index_f"
1461
1462            !    declare optional arguments
1463            integer          :: nm1;   integer          :: dt1
1464            integer,optional:: nm2;   integer,optional:: dt2
1465            integer,optional:: nm3;   integer,optional:: dt3
1466            integer,optional:: nm4;   integer,optional:: dt4
1467            integer,optional:: nm5;   integer,optional:: dt5
1468            integer,optional:: nm6;   integer,optional:: dt6
1469            integer,optional:: nm7;   integer,optional:: dt7
1470            integer,optional:: nm8;   integer,optional:: dt8
1471            integer,optional:: nm9;   integer,optional:: dt9
1472            integer,optional:: nm10;  integer,optional:: dt10
1473            integer,optional:: nm11;  integer,optional:: dt11
1474            integer,optional:: nm12;  integer,optional:: dt12
1475            integer,optional:: nm13;  integer,optional:: dt13
1476            integer,optional:: nm14;  integer,optional:: dt14
1477            integer,optional:: nm15;  integer,optional:: dt15
1478            integer ::k1;
1479            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1480            !
1481
1482            if(k1.ne.0.and.k1.ne.1) then
1483               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1484               if(errcheck_mltpl(error,msg,1,"", nm1 ).lt.KIM_STATUS_OK) return
1485            endif
1486            if(k1.eq.1) then
1487                dt1 = kim_api_get_compute_by_index_f(kimmdl,nm1,error);
1488                if (errcheck_mltpl(error,msg,1, "", nm1).lt.KIM_STATUS_OK) return
1489            endif
1490
1491            !check rest of the arguments
1492            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1493            if( present(nm2 ).and.(.not.present(dt2 ))) then
1494                if(errcheck_mltpl(error,msg,2, "",nm2 ).lt.KIM_STATUS_OK) return
1495            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1496                if(errcheck_mltpl(error,msg,3, "",nm3 ).lt.KIM_STATUS_OK) return
1497            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1498                if(errcheck_mltpl(error,msg,4, "",nm4 ).lt.KIM_STATUS_OK) return
1499            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1500                if(errcheck_mltpl(error,msg,5, "",nm5 ).lt.KIM_STATUS_OK) return
1501            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1502                if(errcheck_mltpl(error,msg,6, "",nm6 ).lt.KIM_STATUS_OK) return
1503            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1504                if(errcheck_mltpl(error,msg,7, "",nm7 ).lt.KIM_STATUS_OK) return
1505            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1506                if(errcheck_mltpl(error,msg,8, "",nm8 ).lt.KIM_STATUS_OK) return
1507            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1508                if(errcheck_mltpl(error,msg,9, "",nm9 ).lt.KIM_STATUS_OK) return
1509            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1510                if(errcheck_mltpl(error,msg,10, "",nm10 ).lt.KIM_STATUS_OK) return
1511            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1512                if(errcheck_mltpl(error,msg,11, "",nm11 ).lt.KIM_STATUS_OK) return
1513            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1514                if(errcheck_mltpl(error,msg,12, "",nm12 ).lt.KIM_STATUS_OK) return
1515            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1516                if(errcheck_mltpl(error,msg,13, "",nm13 ).lt.KIM_STATUS_OK) return
1517            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1518                if(errcheck_mltpl(error,msg,14, "",nm14 ).lt.KIM_STATUS_OK) return
1519            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1520                if(errcheck_mltpl(error,msg,15, "",nm15 ).lt.KIM_STATUS_OK) return
1521            endif
1522
1523            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1524            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1525                if(errcheck_mltpl(error,msg,2,"", nm2 ).lt.KIM_STATUS_OK) return
1526            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1527                if(errcheck_mltpl(error,msg,3,"", nm3 ).lt.KIM_STATUS_OK) return
1528            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1529                if(errcheck_mltpl(error,msg,4,"", nm4 ).lt.KIM_STATUS_OK) return
1530            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1531                if(errcheck_mltpl(error,msg,5,"", nm5 ).lt.KIM_STATUS_OK) return
1532            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1533                if(errcheck_mltpl(error,msg,6,"", nm6 ).lt.KIM_STATUS_OK) return
1534            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1535                if(errcheck_mltpl(error,msg,7,"", nm7 ).lt.KIM_STATUS_OK) return
1536            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1537                if(errcheck_mltpl(error,msg,8,"", nm8 ).lt.KIM_STATUS_OK) return
1538            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1539                if(errcheck_mltpl(error,msg,9,"", nm9 ).lt.KIM_STATUS_OK) return
1540            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1541                if(errcheck_mltpl(error,msg,10,"", nm10 ).lt.KIM_STATUS_OK) return
1542            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1543                if(errcheck_mltpl(error,msg,11,"", nm11 ).lt.KIM_STATUS_OK) return
1544            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1545                if(errcheck_mltpl(error,msg,12,"", nm12 ).lt.KIM_STATUS_OK) return
1546            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1547                if(errcheck_mltpl(error,msg,13,"", nm13 ).lt.KIM_STATUS_OK) return
1548            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1549                if(errcheck_mltpl(error,msg,14,"", nm14 ).lt.KIM_STATUS_OK) return
1550            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1551                if(errcheck_mltpl(error,msg,15,"", nm15 ).lt.KIM_STATUS_OK) return
1552            endif
1553
1554            !process arguments
1555            error =KIM_STATUS_OK
1556 if(present(nm2).and.k2.eq.1) dt2= kim_api_get_compute_by_index_f(kimmdl,nm2,error);
1557            if(errcheck_mltpl(error,msg,2,"",nm2).lt.KIM_STATUS_OK) return
1558 if(present(nm3).and.k3.eq.1) dt3= kim_api_get_compute_by_index_f(kimmdl,nm3,error);
1559            if(errcheck_mltpl(error,msg,3,"",nm3).lt.KIM_STATUS_OK) return
1560 if(present(nm4).and.k4.eq.1) dt4= kim_api_get_compute_by_index_f(kimmdl,nm4,error);
1561            if(errcheck_mltpl(error,msg,4,"",nm4).lt.KIM_STATUS_OK) return
1562 if(present(nm5).and.k5.eq.1) dt5= kim_api_get_compute_by_index_f(kimmdl,nm5,error);
1563            if(errcheck_mltpl(error,msg,5,"",nm5).lt.KIM_STATUS_OK) return
1564 if(present(nm6).and.k6.eq.1) dt6= kim_api_get_compute_by_index_f(kimmdl,nm6,error);
1565            if(errcheck_mltpl(error,msg,6,"",nm6).lt.KIM_STATUS_OK) return
1566 if(present(nm7).and.k7.eq.1) dt7= kim_api_get_compute_by_index_f(kimmdl,nm7,error);
1567            if(errcheck_mltpl(error,msg,7,"",nm7).lt.KIM_STATUS_OK) return
1568 if(present(nm8).and.k8.eq.1) dt8= kim_api_get_compute_by_index_f(kimmdl,nm8,error);
1569            if(errcheck_mltpl(error,msg,8,"",nm8).lt.KIM_STATUS_OK) return
1570 if(present(nm9).and.k9.eq.1) dt9= kim_api_get_compute_by_index_f(kimmdl,nm9,error);
1571            if(errcheck_mltpl(error,msg,9,"",nm9).lt.KIM_STATUS_OK) return
1572 if(present(nm10).and.k10.eq.1) dt10= kim_api_get_compute_by_index_f(kimmdl,nm10,error);
1573            if(errcheck_mltpl(error,msg,10,"",nm10).lt.KIM_STATUS_OK) return
1574 if(present(nm11).and.k11.eq.1) dt11= kim_api_get_compute_by_index_f(kimmdl,nm11,error);
1575            if(errcheck_mltpl(error,msg,11,"",nm11).lt.KIM_STATUS_OK) return
1576 if(present(nm12).and.k12.eq.1) dt12= kim_api_get_compute_by_index_f(kimmdl,nm12,error);
1577            if(errcheck_mltpl(error,msg,12,"",nm12).lt.KIM_STATUS_OK) return
1578 if(present(nm13).and.k13.eq.1) dt13= kim_api_get_compute_by_index_f(kimmdl,nm13,error);
1579            if(errcheck_mltpl(error,msg,13,"",nm13).lt.KIM_STATUS_OK) return
1580 if(present(nm14).and.k14.eq.1) dt14= kim_api_get_compute_by_index_f(kimmdl,nm14,error);
1581            if(errcheck_mltpl(error,msg,14,"",nm14).lt.KIM_STATUS_OK) return
1582 if(present(nm15).and.k15.eq.1) dt15= kim_api_get_compute_by_index_f(kimmdl,nm15,error);
1583            if(errcheck_mltpl(error,msg,15,"",nm15).lt.KIM_STATUS_OK) return
1584        end subroutine kim_api_getm_compute_by_index_f
1585
1586        subroutine kim_api_getm_index_f(kimmdl,error,  nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1587         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1588         nm14,dt14,k14, nm15,dt15,k15)
1589            integer(kind=kim_intptr) :: kimmdl;  integer error
1590
1591            character(len=40) ::msg="kim_api_getm_index_f"
1592
1593            !    declare optional arguments
1594            character(len=*)          :: nm1;   integer          :: dt1
1595            character(len=*),optional:: nm2;   integer,optional:: dt2
1596            character(len=*),optional:: nm3;   integer,optional:: dt3
1597            character(len=*),optional:: nm4;   integer,optional:: dt4
1598            character(len=*),optional:: nm5;   integer,optional:: dt5
1599            character(len=*),optional:: nm6;   integer,optional:: dt6
1600            character(len=*),optional:: nm7;   integer,optional:: dt7
1601            character(len=*),optional:: nm8;   integer,optional:: dt8
1602            character(len=*),optional:: nm9;   integer,optional:: dt9
1603            character(len=*),optional:: nm10;  integer,optional:: dt10
1604            character(len=*),optional:: nm11;  integer,optional:: dt11
1605            character(len=*),optional:: nm12;  integer,optional:: dt12
1606            character(len=*),optional:: nm13;  integer,optional:: dt13
1607            character(len=*),optional:: nm14;  integer,optional:: dt14
1608            character(len=*),optional:: nm15;  integer,optional:: dt15
1609            integer ::k1;
1610            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1611            !
1612            if(k1.ne.0.and.k1.ne.1) then
1613               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1614               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1, nm1 ).lt.KIM_STATUS_OK) return
1615            endif
1616            if(k1.eq.1) then
1617                dt1 = kim_api_get_index_f(kimmdl,nm1,error);
1618                if (errcheck_mltpl(error,msg,1, nm1).lt.KIM_STATUS_OK) return
1619            endif
1620
1621            !check rest of the arguments
1622            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1623            if( present(nm2 ).and.(.not.present(dt2 ))) then
1624                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1625            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1626                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1627            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1628                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1629            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1630                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1631            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1632                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1633            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1634                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1635            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1636                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1637            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1638                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1639            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1640                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1641            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1642                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1643            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1644                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1645            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1646                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1647            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1648                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1649            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1650                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1651            endif
1652
1653            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1654            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1655                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1656            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1657                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1658            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1659                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1660            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1661                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1662            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1663                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1664            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1665                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1666            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1667                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1668            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1669                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1670            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1671                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1672            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1673                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1674            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1675                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1676            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1677                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1678            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1679                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1680            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1681                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1682            endif
1683
1684            !process arguments
1685            error = KIM_STATUS_OK
1686 if(present(nm2).and.k2.eq.1) dt2= kim_api_get_index_f(kimmdl,nm2,error);
1687            if(errcheck_mltpl(error,msg,2,nm2).lt.KIM_STATUS_OK) return
1688 if(present(nm3).and.k3.eq.1) dt3= kim_api_get_index_f(kimmdl,nm3,error);
1689            if(errcheck_mltpl(error,msg,3,nm3).lt.KIM_STATUS_OK) return
1690 if(present(nm4).and.k4.eq.1) dt4= kim_api_get_index_f(kimmdl,nm4,error);
1691            if(errcheck_mltpl(error,msg,4,nm4).lt.KIM_STATUS_OK) return
1692 if(present(nm5).and.k5.eq.1) dt5= kim_api_get_index_f(kimmdl,nm5,error);
1693            if(errcheck_mltpl(error,msg,5,nm5).lt.KIM_STATUS_OK) return
1694 if(present(nm6).and.k6.eq.1) dt6= kim_api_get_index_f(kimmdl,nm6,error);
1695            if(errcheck_mltpl(error,msg,6,nm6).lt.KIM_STATUS_OK) return
1696 if(present(nm7).and.k7.eq.1) dt7= kim_api_get_index_f(kimmdl,nm7,error);
1697            if(errcheck_mltpl(error,msg,7,nm7).lt.KIM_STATUS_OK) return
1698 if(present(nm8).and.k8.eq.1) dt8= kim_api_get_index_f(kimmdl,nm8,error);
1699            if(errcheck_mltpl(error,msg,8,nm8).lt.KIM_STATUS_OK) return
1700 if(present(nm9).and.k9.eq.1) dt9= kim_api_get_index_f(kimmdl,nm9,error);
1701            if(errcheck_mltpl(error,msg,9,nm9).lt.KIM_STATUS_OK) return
1702 if(present(nm10).and.k10.eq.1) dt10= kim_api_get_index_f(kimmdl,nm10,error);
1703            if(errcheck_mltpl(error,msg,10,nm10).lt.KIM_STATUS_OK) return
1704 if(present(nm11).and.k11.eq.1) dt11= kim_api_get_index_f(kimmdl,nm11,error);
1705            if(errcheck_mltpl(error,msg,11,nm11).lt.KIM_STATUS_OK) return
1706 if(present(nm12).and.k12.eq.1) dt12= kim_api_get_index_f(kimmdl,nm12,error);
1707            if(errcheck_mltpl(error,msg,12,nm12).lt.KIM_STATUS_OK) return
1708 if(present(nm13).and.k13.eq.1) dt13= kim_api_get_index_f(kimmdl,nm13,error);
1709            if(errcheck_mltpl(error,msg,13,nm13).lt.KIM_STATUS_OK) return
1710 if(present(nm14).and.k14.eq.1) dt14= kim_api_get_index_f(kimmdl,nm14,error);
1711            if(errcheck_mltpl(error,msg,14,nm14).lt.KIM_STATUS_OK) return
1712 if(present(nm15).and.k15.eq.1) dt15= kim_api_get_index_f(kimmdl,nm15,error);
1713            if(errcheck_mltpl(error,msg,15,nm15).lt.KIM_STATUS_OK) return
1714        end subroutine kim_api_getm_index_f
1715
1716         subroutine kim_api_setm_compute_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1717         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1718         nm14,dt14,k14, nm15,dt15,k15)
1719            integer(kind=kim_intptr) :: kimmdl;  integer error
1720
1721            character(len=40) ::msg="kim_api_setm_compute_f"
1722
1723            !    declare optional arguments
1724            character(len=*)          :: nm1;   integer          :: dt1
1725            character(len=*),optional:: nm2;   integer,optional:: dt2
1726            character(len=*),optional:: nm3;   integer,optional:: dt3
1727            character(len=*),optional:: nm4;   integer,optional:: dt4
1728            character(len=*),optional:: nm5;   integer,optional:: dt5
1729            character(len=*),optional:: nm6;   integer,optional:: dt6
1730            character(len=*),optional:: nm7;   integer,optional:: dt7
1731            character(len=*),optional:: nm8;   integer,optional:: dt8
1732            character(len=*),optional:: nm9;   integer,optional:: dt9
1733            character(len=*),optional:: nm10;  integer,optional:: dt10
1734            character(len=*),optional:: nm11;  integer,optional:: dt11
1735            character(len=*),optional:: nm12;  integer,optional:: dt12
1736            character(len=*),optional:: nm13;  integer,optional:: dt13
1737            character(len=*),optional:: nm14;  integer,optional:: dt14
1738            character(len=*),optional:: nm15;  integer,optional:: dt15
1739            integer ::k1;
1740            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1741            !
1742            if(k1.ne.0.and.k1.ne.1) then
1743                error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1744               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1, nm1 ).lt.KIM_STATUS_OK) return
1745            endif
1746            if(k1.eq.1) then
1747                call kim_api_set_compute_f(kimmdl,nm1,dt1,error);
1748                if (errcheck_mltpl(error,msg,1, nm1).lt.KIM_STATUS_OK) return
1749            endif
1750
1751            !check rest of the arguments
1752            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1753            if( present(nm2 ).and.(.not.present(dt2 ))) then
1754                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1755            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1756                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1757            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1758                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1759            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1760                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1761            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1762                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1763            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1764                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1765            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1766                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1767            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1768                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1769            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1770                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1771            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1772                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1773            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1774                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1775            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1776                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1777            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1778                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1779            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1780                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1781            endif
1782
1783            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1784            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1785                if(errcheck_mltpl(error,msg,2, nm2 ).lt.KIM_STATUS_OK) return
1786            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1787                if(errcheck_mltpl(error,msg,3, nm3 ).lt.KIM_STATUS_OK) return
1788            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1789                if(errcheck_mltpl(error,msg,4, nm4 ).lt.KIM_STATUS_OK) return
1790            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1791                if(errcheck_mltpl(error,msg,5, nm5 ).lt.KIM_STATUS_OK) return
1792            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1793                if(errcheck_mltpl(error,msg,6, nm6 ).lt.KIM_STATUS_OK) return
1794            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1795                if(errcheck_mltpl(error,msg,7, nm7 ).lt.KIM_STATUS_OK) return
1796            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1797                if(errcheck_mltpl(error,msg,8, nm8 ).lt.KIM_STATUS_OK) return
1798            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1799                if(errcheck_mltpl(error,msg,9, nm9 ).lt.KIM_STATUS_OK) return
1800            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1801                if(errcheck_mltpl(error,msg,10, nm10 ).lt.KIM_STATUS_OK) return
1802            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1803                if(errcheck_mltpl(error,msg,11, nm11 ).lt.KIM_STATUS_OK) return
1804            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1805                if(errcheck_mltpl(error,msg,12, nm12 ).lt.KIM_STATUS_OK) return
1806            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1807                if(errcheck_mltpl(error,msg,13, nm13 ).lt.KIM_STATUS_OK) return
1808            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1809                if(errcheck_mltpl(error,msg,14, nm14 ).lt.KIM_STATUS_OK) return
1810            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1811                if(errcheck_mltpl(error,msg,15, nm15 ).lt.KIM_STATUS_OK) return
1812            endif
1813
1814            !process arguments
1815            error = KIM_STATUS_OK
1816 if(present(nm2).and.k2.eq.1) call kim_api_set_compute_f(kimmdl,nm2,dt2,error);
1817            if(errcheck_mltpl(error,msg,2,nm2).lt.KIM_STATUS_OK) return
1818 if(present(nm3).and.k3.eq.1) call kim_api_set_compute_f(kimmdl,nm3,dt3,error);
1819            if(errcheck_mltpl(error,msg,3,nm3).lt.KIM_STATUS_OK) return
1820 if(present(nm4).and.k4.eq.1) call kim_api_set_compute_f(kimmdl,nm4,dt4,error);
1821            if(errcheck_mltpl(error,msg,4,nm4).lt.KIM_STATUS_OK) return
1822 if(present(nm5).and.k5.eq.1) call kim_api_set_compute_f(kimmdl,nm5,dt5,error);
1823            if(errcheck_mltpl(error,msg,5,nm5).lt.KIM_STATUS_OK) return
1824 if(present(nm6).and.k6.eq.1) call kim_api_set_compute_f(kimmdl,nm6,dt6,error);
1825            if(errcheck_mltpl(error,msg,6,nm6).lt.KIM_STATUS_OK) return
1826 if(present(nm7).and.k7.eq.1) call kim_api_set_compute_f(kimmdl,nm7,dt7,error);
1827            if(errcheck_mltpl(error,msg,7,nm7).lt.KIM_STATUS_OK) return
1828 if(present(nm8).and.k8.eq.1) call kim_api_set_compute_f(kimmdl,nm8,dt8,error);
1829            if(errcheck_mltpl(error,msg,8,nm8).lt.KIM_STATUS_OK) return
1830 if(present(nm9).and.k9.eq.1) call kim_api_set_compute_f(kimmdl,nm9,dt9,error);
1831            if(errcheck_mltpl(error,msg,9,nm9).lt.KIM_STATUS_OK) return
1832 if(present(nm10).and.k10.eq.1) call kim_api_set_compute_f(kimmdl,nm10,dt10,error);
1833            if(errcheck_mltpl(error,msg,10,nm10).lt.KIM_STATUS_OK) return
1834 if(present(nm11).and.k11.eq.1) call kim_api_set_compute_f(kimmdl,nm11,dt11,error);
1835            if(errcheck_mltpl(error,msg,11,nm11).lt.KIM_STATUS_OK) return
1836 if(present(nm12).and.k12.eq.1) call kim_api_set_compute_f(kimmdl,nm12,dt12,error);
1837            if(errcheck_mltpl(error,msg,12,nm12).lt.KIM_STATUS_OK) return
1838 if(present(nm13).and.k13.eq.1) call kim_api_set_compute_f(kimmdl,nm13,dt13,error);
1839            if(errcheck_mltpl(error,msg,13,nm13).lt.KIM_STATUS_OK) return
1840 if(present(nm14).and.k14.eq.1) call kim_api_set_compute_f(kimmdl,nm14,dt14,error);
1841            if(errcheck_mltpl(error,msg,14,nm14).lt.KIM_STATUS_OK) return
1842 if(present(nm15).and.k15.eq.1) call kim_api_set_compute_f(kimmdl,nm15,dt15,error);
1843            if(errcheck_mltpl(error,msg,15,nm15).lt.KIM_STATUS_OK) return
1844         end subroutine kim_api_setm_compute_f
1845
1846        subroutine kim_api_setm_compute_by_index_f(kimmdl,error, nm1,dt1,k1, nm2,dt2,k2, nm3,dt3,k3, nm4,dt4,k4, nm5,dt5,k5,&
1847         nm6,dt6,k6, nm7,dt7,k7, nm8,dt8,k8, nm9,dt9,k9, nm10,dt10,k10, nm11,dt11,k11, nm12,dt12,k12, nm13,dt13,k13,&
1848         nm14,dt14,k14, nm15,dt15,k15)
1849            integer(kind=kim_intptr) :: kimmdl;  integer error
1850
1851            character(len=40) ::msg="kim_api_set_compute_byi_mltpl_f"
1852
1853            !    declare optional arguments
1854            integer          :: nm1;   integer          :: dt1
1855            integer,optional:: nm2;   integer,optional:: dt2
1856            integer,optional:: nm3;   integer,optional:: dt3
1857            integer,optional:: nm4;   integer,optional:: dt4
1858            integer,optional:: nm5;   integer,optional:: dt5
1859            integer,optional:: nm6;   integer,optional:: dt6
1860            integer,optional:: nm7;   integer,optional:: dt7
1861            integer,optional:: nm8;   integer,optional:: dt8
1862            integer,optional:: nm9;   integer,optional:: dt9
1863            integer,optional:: nm10;  integer,optional:: dt10
1864            integer,optional:: nm11;  integer,optional:: dt11
1865            integer,optional:: nm12;  integer,optional:: dt12
1866            integer,optional:: nm13;  integer,optional:: dt13
1867            integer,optional:: nm14;  integer,optional:: dt14
1868            integer,optional:: nm15;  integer,optional:: dt15
1869            integer ::k1;
1870            integer,optional::k2,k3,k4,k5,k6,k7,k8,k9,k10,k11,k12,k13,k14,k15
1871            !
1872            if(k1.ne.0.and.k1.ne.1) then
1873               error=KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1874               if(errcheck_mltpl(KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY,msg,1,"", nm1 ).lt.KIM_STATUS_OK) return
1875            endif
1876            if(k1.eq.1) then
1877                call kim_api_set_compute_by_index_f(kimmdl,nm1,dt1,error);
1878                if (errcheck_mltpl(error,msg,1, "", nm1).lt.KIM_STATUS_OK) return
1879            endif
1880
1881            !check rest of the arguments
1882            error = KIM_STATUS_WRONG_MULTIPLE_ARGS
1883            if( present(nm2 ).and.(.not.present(dt2 ))) then
1884                if(errcheck_mltpl(error,msg,2, "",nm2 ).lt.KIM_STATUS_OK) return
1885            elseif( present(nm3 ).and.(.not.present(dt3 ))) then
1886                if(errcheck_mltpl(error,msg,3, "",nm3 ).lt.KIM_STATUS_OK) return
1887            elseif( present(nm4 ).and.(.not.present(dt4 ))) then
1888                if(errcheck_mltpl(error,msg,4, "",nm4 ).lt.KIM_STATUS_OK) return
1889            elseif( present(nm5 ).and.(.not.present(dt5 ))) then
1890                if(errcheck_mltpl(error,msg,5, "",nm5 ).lt.KIM_STATUS_OK) return
1891            elseif( present(nm6 ).and.(.not.present(dt6 ))) then
1892                if(errcheck_mltpl(error,msg,6, "",nm6 ).lt.KIM_STATUS_OK) return
1893            elseif( present(nm7 ).and.(.not.present(dt7 ))) then
1894                if(errcheck_mltpl(error,msg,7, "",nm7 ).lt.KIM_STATUS_OK) return
1895            elseif( present(nm8 ).and.(.not.present(dt8 ))) then
1896                if(errcheck_mltpl(error,msg,8, "",nm8 ).lt.KIM_STATUS_OK) return
1897            elseif( present(nm9 ).and.(.not.present(dt9 ))) then
1898                if(errcheck_mltpl(error,msg,9, "",nm9 ).lt.KIM_STATUS_OK) return
1899            elseif( present(nm10 ).and.(.not.present(dt10 ))) then
1900                if(errcheck_mltpl(error,msg,10, "",nm10 ).lt.KIM_STATUS_OK) return
1901            elseif( present(nm11 ).and.(.not.present(dt11 ))) then
1902                if(errcheck_mltpl(error,msg,11, "",nm11 ).lt.KIM_STATUS_OK) return
1903            elseif( present(nm12 ).and.(.not.present(dt12 ))) then
1904                if(errcheck_mltpl(error,msg,12, "",nm12 ).lt.KIM_STATUS_OK) return
1905            elseif( present(nm13 ).and.(.not.present(dt13 ))) then
1906                if(errcheck_mltpl(error,msg,13, "",nm13 ).lt.KIM_STATUS_OK) return
1907            elseif( present(nm14 ).and.(.not.present(dt14 ))) then
1908                if(errcheck_mltpl(error,msg,14, "",nm14 ).lt.KIM_STATUS_OK) return
1909            elseif( present(nm15 ).and.(.not.present(dt15 ))) then
1910                if(errcheck_mltpl(error,msg,15, "",nm15 ).lt.KIM_STATUS_OK) return
1911            endif
1912
1913            error = KIM_STATUS_WRONG_GROUP_ARGUMENT_KEY
1914            if(present(k2).and.(k2.ne.0.and.k2.ne.1))then
1915                if(errcheck_mltpl(error,msg,2,"", nm2 ).lt.KIM_STATUS_OK) return
1916            elseif(present(k3).and.(k3.ne.0.and.k3.ne.1))then
1917                if(errcheck_mltpl(error,msg,3,"", nm3 ).lt.KIM_STATUS_OK) return
1918            elseif(present(k4).and.(k4.ne.0.and.k4.ne.1))then
1919                if(errcheck_mltpl(error,msg,4,"", nm4 ).lt.KIM_STATUS_OK) return
1920            elseif(present(k5).and.(k5.ne.0.and.k5.ne.1))then
1921                if(errcheck_mltpl(error,msg,5,"", nm5 ).lt.KIM_STATUS_OK) return
1922            elseif(present(k6).and.(k6.ne.0.and.k6.ne.1))then
1923                if(errcheck_mltpl(error,msg,6,"", nm6 ).lt.KIM_STATUS_OK) return
1924            elseif(present(k7).and.(k7.ne.0.and.k7.ne.1))then
1925                if(errcheck_mltpl(error,msg,7,"", nm7 ).lt.KIM_STATUS_OK) return
1926            elseif(present(k8).and.(k8.ne.0.and.k8.ne.1))then
1927                if(errcheck_mltpl(error,msg,8,"", nm8 ).lt.KIM_STATUS_OK) return
1928            elseif(present(k9).and.(k9.ne.0.and.k9.ne.1))then
1929                if(errcheck_mltpl(error,msg,9,"", nm9 ).lt.KIM_STATUS_OK) return
1930            elseif(present(k10).and.(k10.ne.0.and.k10.ne.1))then
1931                if(errcheck_mltpl(error,msg,10,"", nm10 ).lt.KIM_STATUS_OK) return
1932            elseif(present(k11).and.(k11.ne.0.and.k11.ne.1))then
1933                if(errcheck_mltpl(error,msg,11,"", nm11 ).lt.KIM_STATUS_OK) return
1934            elseif(present(k12).and.(k12.ne.0.and.k12.ne.1))then
1935                if(errcheck_mltpl(error,msg,12,"", nm12 ).lt.KIM_STATUS_OK) return
1936            elseif(present(k13).and.(k13.ne.0.and.k13.ne.1))then
1937                if(errcheck_mltpl(error,msg,13,"", nm13 ).lt.KIM_STATUS_OK) return
1938            elseif(present(k14).and.(k14.ne.0.and.k14.ne.1))then
1939                if(errcheck_mltpl(error,msg,14,"", nm14 ).lt.KIM_STATUS_OK) return
1940            elseif(present(k15).and.(k15.ne.0.and.k15.ne.1))then
1941                if(errcheck_mltpl(error,msg,15,"", nm15 ).lt.KIM_STATUS_OK) return
1942            endif
1943
1944            !process arguments
1945            error = KIM_STATUS_OK
1946 if(present(nm2).and.k2.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm2,dt2,error);
1947            if(errcheck_mltpl(error,msg,2,"",nm2).lt.KIM_STATUS_OK) return
1948 if(present(nm3).and.k3.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm3,dt3,error);
1949            if(errcheck_mltpl(error,msg,3,"",nm3).lt.KIM_STATUS_OK) return
1950 if(present(nm4).and.k4.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm4,dt4,error);
1951            if(errcheck_mltpl(error,msg,4,"",nm4).lt.KIM_STATUS_OK) return
1952 if(present(nm5).and.k5.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm5,dt5,error);
1953            if(errcheck_mltpl(error,msg,5,"",nm5).lt.KIM_STATUS_OK) return
1954 if(present(nm6).and.k6.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm6,dt6,error);
1955            if(errcheck_mltpl(error,msg,6,"",nm6).lt.KIM_STATUS_OK) return
1956 if(present(nm7).and.k7.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm7,dt7,error);
1957            if(errcheck_mltpl(error,msg,7,"",nm7).lt.KIM_STATUS_OK) return
1958 if(present(nm8).and.k8.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm8,dt8,error);
1959            if(errcheck_mltpl(error,msg,8,"",nm8).lt.KIM_STATUS_OK) return
1960 if(present(nm9).and.k9.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm9,dt9,error);
1961            if(errcheck_mltpl(error,msg,9,"",nm9).lt.KIM_STATUS_OK) return
1962 if(present(nm10).and.k10.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm10,dt10,error);
1963            if(errcheck_mltpl(error,msg,10,"",nm10).lt.KIM_STATUS_OK) return
1964 if(present(nm11).and.k11.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm11,dt11,error);
1965            if(errcheck_mltpl(error,msg,11,"",nm11).lt.KIM_STATUS_OK) return
1966 if(present(nm12).and.k12.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm12,dt12,error);
1967            if(errcheck_mltpl(error,msg,12,"",nm12).lt.KIM_STATUS_OK) return
1968 if(present(nm13).and.k13.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm13,dt13,error);
1969            if(errcheck_mltpl(error,msg,13,"",nm13).lt.KIM_STATUS_OK) return
1970 if(present(nm14).and.k14.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm14,dt14,error);
1971            if(errcheck_mltpl(error,msg,14,"",nm14).lt.KIM_STATUS_OK) return
1972 if(present(nm15).and.k15.eq.1) call kim_api_set_compute_by_index_f(kimmdl,nm15,dt15,error);
1973            if(errcheck_mltpl(error,msg,15,"",nm15).lt.KIM_STATUS_OK) return
1974        end subroutine kim_api_setm_compute_by_index_f
1975
1976
1977       integer function errcheck_mltpl(error,msgfrom,grarg,nm,ind)
1978            integer ::error,grarg; character(len=*)::msgfrom
1979            character(len=*), optional::nm
1980            integer,optional ::ind
1981            errcheck_mltpl = KIM_STATUS_OK
1982            if(error.ge.KIM_STATUS_OK) return
1983            errcheck_mltpl = KIM_STATUS_FAIL
1984            if(present(nm).and.present(ind))then
1985                print*,"failed:", msgfrom,", for group argument ", grarg, " and kim_name ",nm,", kim_index", ind
1986            else if(present(nm)) then
1987                print*,"failed:", msgfrom,", for group argument ", grarg, " and kim_name ",nm
1988            else if(present(ind)) then
1989                print*,"failed:", msgfrom,", for group argument ", grarg, " and kim_index",ind
1990            else
1991                print*,"failed:", msgfrom,", for group argument ", grarg
1992            endif
1993            errcheck_mltpl = KIM_STATUS_FAIL
1994       end function errcheck_mltpl
1995
1996end module kim_api
1997