/[MITgcm]/MITgcm/tools/OAD_support/OAD_active.F90
ViewVC logotype

Annotation of /MITgcm/tools/OAD_support/OAD_active.F90

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download)
Thu Sep 20 23:12:47 2012 UTC (11 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64o, checkpoint64a, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint64n, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64c, checkpoint64g, checkpoint64f, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64, checkpoint65, checkpoint64j, checkpoint64m, checkpoint64l, HEAD
* Merge OAD_support from MITgcm_contrib/heimbach/OpenAD/OAD_support/
  to tools/OAD_support/
* Adjust genmake2 to reflect path change (attempt with ${OADTOOLS})
* Adjust insertTemplateDir.bash to reflect path change
Seems to work.

1 heimbach 1.1 !#########################################################
2     ! This file is part of OpenAD released under the LGPL. #
3     ! The full COPYRIGHT notice can be found in the top #
4     ! level directory of the OpenAD distribution #
5     !#########################################################
6     module OAD_active
7     use w2f__types
8     implicit none
9     private :: runTimeErrorStop, shapeChange
10     public :: active
11     #ifndef TRACE
12     public :: saxpy, sax, zero_deriv, setderiv
13     public :: set_neg_deriv, inc_deriv, dec_deriv
14     #endif
15     public :: oad_convert
16     public :: oad_allocateMatching, oad_allocateShape, oad_shapeTest
17     #ifndef TRACE
18     integer :: count_mult = 0
19     integer :: count_add = 0
20     #endif
21     integer, parameter :: shapeChange=0
22     #ifdef VECTOR
23     integer :: max_deriv_vec_len
24     parameter ( max_deriv_vec_len = 100 )
25     integer :: oad_ctmp_i1, oad_ctmp_i2, oad_ctmp_i3, &
26     oad_ctmp_i4, oad_ctmp_i5, oad_ctmp_i6, oad_ctmp_i7
27     # define VECTOR_DIM , dimension(max_deriv_vec_len)
28     # define VECTOR_LOOP_VAR integer :: i
29     # define VECTOR_LOOP_BEGIN do i=1,max_deriv_vec_len
30     # define VECTOR_LOOP_END end do
31     # define DELEM d(i)
32     #else
33     # define VECTOR_DIM
34     # define VECTOR_LOOP_VAR
35     # define VECTOR_LOOP_BEGIN
36     # define VECTOR_LOOP_END
37     # define DELEM d
38     #endif
39     #ifdef SCALARNDI
40     # define DINIT
41     #else
42     # define DINIT =0.0d0
43     #endif
44    
45     !
46     ! active needs to be a sequence type
47     ! with no initialization
48     !
49     type active
50     sequence
51     real(w2f__8) :: v
52     #ifndef TRACE
53     ! initialization does not work for active variables
54     ! inside of common block, such as in boxmodel
55     ! initialization is required for correct adjoint
56     real(w2f__8) VECTOR_DIM :: d DINIT
57     #endif
58     end type
59     #ifndef TRACE
60     interface saxpy
61     module procedure saxpy_d0_a0_a0
62     module procedure saxpy_l0_a0_a0
63     module procedure saxpy_i0_a0_a0
64     module procedure saxpy_d0_a0_a1
65     module procedure saxpy_d0_a1_a1
66     module procedure saxpy_d0_a2_a2
67     module procedure saxpy_d1_a0_a1
68     module procedure saxpy_d1_a1_a1
69     module procedure saxpy_l1_a1_a1
70     module procedure saxpy_i1_a1_a1
71     module procedure saxpy_a1_a1_a1
72     module procedure saxpy_d2_a0_a2
73     module procedure saxpy_d2_a2_a2
74     # ifndef DEFAULT_R8
75     module procedure saxpy_r0_a0_a0
76     module procedure saxpy_r0_a1_a1
77     module procedure saxpy_r1_a0_a1
78     module procedure saxpy_r1_a1_a1
79     # endif
80     end interface
81    
82     interface setderiv
83     module procedure setderiv_a0_a0
84     module procedure setderiv_a1_a0
85     module procedure setderiv_a1_a1
86     module procedure setderiv_a2_a0
87     module procedure setderiv_a2_a2
88     module procedure setderiv_a3_a3
89     end interface
90    
91     interface set_neg_deriv
92     module procedure set_neg_deriv_a0_a0
93     module procedure set_neg_deriv_a1_a1
94     end interface set_neg_deriv
95    
96     interface inc_deriv
97     module procedure inc_deriv_a0_a0
98     module procedure inc_deriv_a1_a1
99     module procedure inc_deriv_a2_a2
100     end interface inc_deriv
101    
102     interface dec_deriv
103     module procedure dec_deriv_a0_a0
104     module procedure dec_deriv_a1_a1
105     module procedure dec_deriv_a2_a2
106     end interface dec_deriv
107    
108     interface zero_deriv
109     module procedure zero_deriv_a0
110     module procedure zero_deriv_a1
111     module procedure zero_deriv_a2
112     module procedure zero_deriv_a3
113     module procedure zero_deriv_a4
114     end interface
115    
116     interface sax
117     module procedure sax_d0_a0_a0
118     module procedure sax_l0_a0_a0
119     module procedure sax_i0_a0_a0
120     module procedure sax_d0_a0_a1
121     module procedure sax_d0_a0_a2
122     module procedure sax_d0_a1_a1
123     module procedure sax_d0_a2_a2
124     module procedure sax_d0_a3_a3
125     module procedure sax_d1_a0_a1
126     module procedure sax_d1_a1_a1
127     module procedure sax_l1_a1_a1
128     module procedure sax_i1_a1_a1
129     module procedure sax_d2_a0_a2
130     module procedure sax_d2_a2_a2
131     # ifndef DEFAULT_R8
132     module procedure sax_r0_a0_a0
133     # endif
134     end interface
135    
136     #endif
137     interface oad_convert
138     module procedure convert_d0_a0
139     module procedure convert_d1_a1
140     module procedure convert_d2_a2
141     module procedure convert_d3_a3
142     module procedure convert_d4_a4
143     module procedure convert_d5_a5
144     module procedure convert_d6_a6
145     module procedure convert_d7_a7
146     module procedure convert_a0_d0
147     module procedure convert_a1_d1
148     module procedure convert_a2_d2
149     module procedure convert_a3_d3
150     module procedure convert_a4_d4
151     module procedure convert_a5_d5
152     module procedure convert_a6_d6
153     module procedure convert_a7_d7
154     #ifndef DEFAULT_R8
155     module procedure convert_r0_a0
156     module procedure convert_r1_a1
157     module procedure convert_r2_a2
158     module procedure convert_r3_a3
159     module procedure convert_r4_a4
160     module procedure convert_r5_a5
161     module procedure convert_r6_a6
162     module procedure convert_r7_a7
163     module procedure convert_a0_r0
164     module procedure convert_a1_r1
165     module procedure convert_a2_r2
166     module procedure convert_a3_r3
167     module procedure convert_a4_r4
168     module procedure convert_a5_r5
169     module procedure convert_a6_r6
170     module procedure convert_a7_r7
171     #endif
172     end interface
173    
174     interface oad_allocateMatching
175     module procedure allocateMatching_d1_d1
176     module procedure allocateMatching_a1_d1
177     module procedure allocateMatching_d1_a1
178     module procedure allocateMatching_a1_a1
179     module procedure allocateMatching_d2_d2
180     module procedure allocateMatching_a2_d2
181     module procedure allocateMatching_d2_a2
182     module procedure allocateMatching_a2_a2
183     module procedure allocateMatching_d3_d3
184     module procedure allocateMatching_a3_d3
185     module procedure allocateMatching_d3_a3
186     module procedure allocateMatching_a3_a3
187     module procedure allocateMatching_a4_a4
188     module procedure allocateMatching_d4_a4
189     module procedure allocateMatching_d5_d5
190     module procedure allocateMatching_a5_d5
191     module procedure allocateMatching_d5_a5
192     module procedure allocateMatching_a5_a5
193     module procedure allocateMatching_d6_d6
194     module procedure allocateMatching_a6_d6
195     module procedure allocateMatching_d6_a6
196     module procedure allocateMatching_a6_a6
197     #ifndef DEFAULT_R8
198     module procedure allocateMatching_r1_r1
199     module procedure allocateMatching_d1_r1
200     module procedure allocateMatching_r1_d1
201     module procedure allocateMatching_a1_r1
202     module procedure allocateMatching_r1_a1
203     module procedure allocateMatching_a2_r2
204     module procedure allocateMatching_r2_a2
205     #endif
206     end interface
207    
208     interface oad_allocateShape
209     module procedure allocateShape_d1
210     module procedure allocateShape_d2
211     end interface
212    
213     interface oad_shapeTest
214     module procedure shapeTest_a1_d1
215     module procedure shapeTest_a1_a1
216     module procedure shapeTest_d1_a1
217     module procedure shapeTest_a2_d2
218     module procedure shapeTest_a2_a2
219     module procedure shapeTest_d2_a2
220     module procedure shapeTest_a3_a3
221     module procedure shapeTest_d3_a3
222     module procedure shapeTest_a4_a4
223     module procedure shapeTest_d4_a4
224     module procedure shapeTest_a5_a5
225     module procedure shapeTest_d5_a5
226     module procedure shapeTest_a5_d5
227     module procedure shapeTest_a6_a6
228     module procedure shapeTest_d6_a6
229     module procedure shapeTest_a6_d6
230     #ifndef DEFAULT_R8
231     module procedure shapeTest_r1_a1
232     module procedure shapeTest_r2_a2
233     module procedure shapeTest_a1_r1
234     module procedure shapeTest_a2_r2
235     #endif
236     end interface
237    
238     interface runTimeErrorStop
239     module procedure runTimeErrorStopI
240     end interface
241    
242     contains
243     #ifndef TRACE
244     !
245     ! chain rule saxpy to be used in forward and reverse modes
246     !
247     subroutine saxpy_d0_a0_a0(a,x,y)
248     real(w2f__8), intent(in) :: a
249     type(active), intent(in) :: x
250     type(active), intent(inout) :: y
251     VECTOR_LOOP_VAR
252     VECTOR_LOOP_BEGIN
253     y%DELEM=y%DELEM + x%DELEM*a
254     VECTOR_LOOP_END
255     end subroutine
256     subroutine saxpy_l0_a0_a0(a,x,y)
257     integer(kind=w2f__i8), intent(in) :: a
258     type(active), intent(in) :: x
259     type(active), intent(inout) :: y
260     VECTOR_LOOP_VAR
261     VECTOR_LOOP_BEGIN
262     y%DELEM=y%DELEM + x%DELEM*a
263     VECTOR_LOOP_END
264     end subroutine
265     subroutine saxpy_i0_a0_a0(a,x,y)
266     integer(kind=w2f__i4), intent(in) :: a
267     type(active), intent(in) :: x
268     type(active), intent(inout) :: y
269     VECTOR_LOOP_VAR
270     VECTOR_LOOP_BEGIN
271     y%DELEM=y%DELEM + x%DELEM*a
272     VECTOR_LOOP_END
273     end subroutine
274     subroutine saxpy_d0_a0_a1(a,x,y)
275     real(w2f__8), intent(in) :: a
276     type(active), intent(in) :: x
277     type(active), dimension(:), intent(inout) :: y
278     VECTOR_LOOP_VAR
279     VECTOR_LOOP_BEGIN
280     y%DELEM=y%DELEM + x%DELEM*a
281     VECTOR_LOOP_END
282     end subroutine
283     subroutine saxpy_d0_a1_a1(a,x,y)
284     real(w2f__8), intent(in) :: a
285     type(active), dimension(:), intent(in) :: x
286     type(active), dimension(:), intent(inout) :: y
287     VECTOR_LOOP_VAR
288     VECTOR_LOOP_BEGIN
289     y%DELEM=y%DELEM + x%DELEM*a
290     VECTOR_LOOP_END
291     end subroutine
292     subroutine saxpy_d0_a2_a2(a,x,y)
293     real(w2f__8), intent(in) :: a
294     type(active), dimension(:,:), intent(in) :: x
295     type(active), dimension(:,:), intent(inout) :: y
296     VECTOR_LOOP_VAR
297     VECTOR_LOOP_BEGIN
298     y%DELEM=y%DELEM + x%DELEM*a
299     VECTOR_LOOP_END
300     end subroutine
301     subroutine saxpy_d1_a0_a1(a,x,y)
302     real(w2f__8), dimension(:), intent(in) :: a
303     type(active), intent(in) :: x
304     type(active), dimension(:), intent(inout) :: y
305     VECTOR_LOOP_VAR
306     VECTOR_LOOP_BEGIN
307     y%DELEM=y%DELEM+x%DELEM*a
308     VECTOR_LOOP_END
309     end subroutine
310     subroutine saxpy_d1_a1_a1(a,x,y)
311     real(w2f__8), dimension(:), intent(in) :: a
312     type(active), dimension(:), intent(in) :: x
313     type(active), dimension(:), intent(inout) :: y
314     VECTOR_LOOP_VAR
315     VECTOR_LOOP_BEGIN
316     y%DELEM=y%DELEM+x%DELEM*a
317     VECTOR_LOOP_END
318     end subroutine
319     subroutine saxpy_l1_a1_a1(a,x,y)
320     integer(kind=w2f__i8), dimension(:), intent(in) :: a
321     type(active), dimension(:), intent(in) :: x
322     type(active), dimension(:), intent(inout) :: y
323     VECTOR_LOOP_VAR
324     VECTOR_LOOP_BEGIN
325     y%DELEM=y%DELEM+x%DELEM*a
326     VECTOR_LOOP_END
327     end subroutine
328     subroutine saxpy_i1_a1_a1(a,x,y)
329     integer(kind=w2f__i4), dimension(:), intent(in) :: a
330     type(active), dimension(:), intent(in) :: x
331     type(active), dimension(:), intent(inout) :: y
332     VECTOR_LOOP_VAR
333     VECTOR_LOOP_BEGIN
334     y%DELEM=y%DELEM+x%DELEM*a
335     VECTOR_LOOP_END
336     end subroutine
337     subroutine saxpy_a1_a1_a1(a,x,y)
338     type(active), dimension(:), intent(in) :: a
339     type(active), dimension(:), intent(in) :: x
340     type(active), dimension(:), intent(inout) :: y
341     VECTOR_LOOP_VAR
342     VECTOR_LOOP_BEGIN
343     y%DELEM=y%DELEM+x%DELEM*a%v
344     VECTOR_LOOP_END
345     end subroutine
346     subroutine saxpy_d2_a0_a2(a,x,y)
347     real(w2f__8), dimension(:,:), intent(in) :: a
348     type(active), intent(in) :: x
349     type(active), dimension(:,:), intent(inout) :: y
350     VECTOR_LOOP_VAR
351     VECTOR_LOOP_BEGIN
352     y%DELEM=y%DELEM + x%DELEM*a
353     VECTOR_LOOP_END
354     end subroutine
355     subroutine saxpy_d2_a2_a2(a,x,y)
356     real(w2f__8), dimension(:,:), intent(in) :: a
357     type(active), dimension(:,:), intent(in) :: x
358     type(active), dimension(:,:), intent(inout) :: y
359     VECTOR_LOOP_VAR
360     VECTOR_LOOP_BEGIN
361     y%DELEM=y%DELEM + x%DELEM*a
362     VECTOR_LOOP_END
363     end subroutine
364     # ifndef DEFAULT_R8
365     subroutine saxpy_r0_a0_a0(a,x,y)
366     real(w2f__4), intent(in) :: a
367     type(active), intent(in) :: x
368     type(active), intent(inout) :: y
369     VECTOR_LOOP_VAR
370     VECTOR_LOOP_BEGIN
371     y%DELEM=y%DELEM + x%DELEM*a
372     VECTOR_LOOP_END
373     end subroutine
374     subroutine saxpy_r0_a1_a1(a,x,y)
375     real(w2f__4), intent(in) :: a
376     type(active), dimension(:), intent(in) :: x
377     type(active), dimension(:), intent(inout) :: y
378     VECTOR_LOOP_VAR
379     VECTOR_LOOP_BEGIN
380     y%DELEM=y%DELEM + x%DELEM*a
381     VECTOR_LOOP_END
382     end subroutine
383     subroutine saxpy_r1_a0_a1(a,x,y)
384     real(w2f__4), dimension(:), intent(in) :: a
385     type(active), intent(in) :: x
386     type(active), dimension(:), intent(inout) :: y
387     VECTOR_LOOP_VAR
388     VECTOR_LOOP_BEGIN
389     y%DELEM=y%DELEM + x%DELEM*a
390     VECTOR_LOOP_END
391     end subroutine
392     subroutine saxpy_r1_a1_a1(a,x,y)
393     real(w2f__4), dimension(:), intent(in) :: a
394     type(active), dimension(:), intent(in) :: x
395     type(active), dimension(:), intent(inout) :: y
396     VECTOR_LOOP_VAR
397     VECTOR_LOOP_BEGIN
398     y%DELEM=y%DELEM+x%DELEM*a
399     VECTOR_LOOP_END
400     end subroutine
401     # endif
402     !
403     ! chain rule saxpy to be used in forward and reverse modes
404     ! derivative component of y is equal to zero initially
405     ! note: y needs to be inout as otherwise value component gets
406     ! zeroed out
407     !
408     subroutine sax_d0_a0_a0(a,x,y)
409     real(w2f__8), intent(in) :: a
410     type(active), intent(in) :: x
411     type(active), intent(inout) :: y
412     VECTOR_LOOP_VAR
413     VECTOR_LOOP_BEGIN
414     y%DELEM=x%DELEM*a
415     VECTOR_LOOP_END
416     end subroutine
417     subroutine sax_l0_a0_a0(a,x,y)
418     integer(kind=w2f__i8), intent(in) :: a
419     type(active), intent(in) :: x
420     type(active), intent(inout) :: y
421     VECTOR_LOOP_VAR
422     VECTOR_LOOP_BEGIN
423     y%DELEM=x%DELEM*a
424     VECTOR_LOOP_END
425     end subroutine
426     subroutine sax_i0_a0_a0(a,x,y)
427     integer(kind=w2f__i4), intent(in) :: a
428     type(active), intent(in) :: x
429     type(active), intent(inout) :: y
430     VECTOR_LOOP_VAR
431     VECTOR_LOOP_BEGIN
432     y%DELEM=x%DELEM*a
433     VECTOR_LOOP_END
434     end subroutine
435     subroutine sax_d0_a0_a1(a,x,y)
436     real(w2f__8), intent(in) :: a
437     type(active), intent(in) :: x
438     type(active), dimension(:), intent(inout) :: y
439     VECTOR_LOOP_VAR
440     VECTOR_LOOP_BEGIN
441     y%DELEM=x%DELEM*a
442     VECTOR_LOOP_END
443     end subroutine
444     subroutine sax_d0_a0_a2(a,x,y)
445     real(w2f__8), intent(in) :: a
446     type(active), intent(in) :: x
447     type(active), dimension(:,:), intent(inout) :: y
448     VECTOR_LOOP_VAR
449     VECTOR_LOOP_BEGIN
450     y%DELEM = x%DELEM*a
451     VECTOR_LOOP_END
452     end subroutine
453     subroutine sax_d0_a1_a1(a,x,y)
454     real(w2f__8), intent(in) :: a
455     type(active), dimension(:), intent(in) :: x
456     type(active), dimension(:), intent(inout) :: y
457     VECTOR_LOOP_VAR
458     VECTOR_LOOP_BEGIN
459     y%DELEM=x%DELEM*a
460     VECTOR_LOOP_END
461     end subroutine
462     subroutine sax_d0_a2_a2(a,x,y)
463     real(w2f__8), intent(in) :: a
464     type(active), dimension(:,:), intent(in) :: x
465     type(active), dimension(:,:), intent(inout) :: y
466     VECTOR_LOOP_VAR
467     VECTOR_LOOP_BEGIN
468     y%DELEM=x%DELEM*a
469     VECTOR_LOOP_END
470     end subroutine
471     subroutine sax_d0_a3_a3(a,x,y)
472     real(w2f__8), intent(in) :: a
473     type(active), dimension(:,:,:), intent(in) :: x
474     type(active), dimension(:,:,:), intent(inout) :: y
475     VECTOR_LOOP_VAR
476     VECTOR_LOOP_BEGIN
477     y%DELEM=x%DELEM*a
478     VECTOR_LOOP_END
479     end subroutine
480     subroutine sax_d1_a0_a1(a,x,y)
481     real(w2f__8), dimension(:), intent(in) :: a
482     type(active), intent(in) :: x
483     type(active), dimension(:), intent(inout) :: y
484     VECTOR_LOOP_VAR
485     VECTOR_LOOP_BEGIN
486     y%DELEM=x%DELEM*a
487     VECTOR_LOOP_END
488     end subroutine
489     subroutine sax_d1_a1_a1(a,x,y)
490     real(w2f__8), dimension(:), intent(in) :: a
491     type(active), dimension(:), intent(in) :: x
492     type(active), dimension(:), intent(inout) :: y
493     VECTOR_LOOP_VAR
494     VECTOR_LOOP_BEGIN
495     y%DELEM=x%DELEM*a
496     VECTOR_LOOP_END
497     end subroutine
498     subroutine sax_l1_a1_a1(a,x,y)
499     integer(kind=w2f__i8), dimension(:), intent(in) :: a
500     type(active), dimension(:), intent(in) :: x
501     type(active), dimension(:), intent(inout) :: y
502     VECTOR_LOOP_VAR
503     VECTOR_LOOP_BEGIN
504     y%DELEM=x%DELEM*a
505     VECTOR_LOOP_END
506     end subroutine
507     subroutine sax_i1_a1_a1(a,x,y)
508     integer(kind=w2f__i4), dimension(:), intent(in) :: a
509     type(active), dimension(:), intent(in) :: x
510     type(active), dimension(:), intent(inout) :: y
511     VECTOR_LOOP_VAR
512     VECTOR_LOOP_BEGIN
513     y%DELEM=x%DELEM*a
514     VECTOR_LOOP_END
515     end subroutine
516     subroutine sax_d2_a0_a2(a,x,y)
517     real(w2f__8), dimension(:,:), intent(in) :: a
518     type(active), intent(in) :: x
519     type(active), dimension(:,:), intent(inout) :: y
520     VECTOR_LOOP_VAR
521     VECTOR_LOOP_BEGIN
522     y%DELEM=x%DELEM*a
523     VECTOR_LOOP_END
524     end subroutine
525     subroutine sax_d2_a2_a2(a,x,y)
526     real(w2f__8), dimension(:,:), intent(in) :: a
527     type(active), dimension(:,:), intent(in) :: x
528     type(active), dimension(:,:), intent(inout) :: y
529     VECTOR_LOOP_VAR
530     VECTOR_LOOP_BEGIN
531     y%DELEM=x%DELEM*a
532     VECTOR_LOOP_END
533     end subroutine
534     # ifndef DEFAULT_R8
535     subroutine sax_r0_a0_a0(a,x,y)
536     real(w2f__4), intent(in) :: a
537     type(active), intent(in) :: x
538     type(active), intent(inout) :: y
539     VECTOR_LOOP_VAR
540     VECTOR_LOOP_BEGIN
541     y%DELEM=x%DELEM*a
542     VECTOR_LOOP_END
543     end subroutine
544     # endif
545     !
546     ! set derivative of y to be equal to derivative of x
547     ! note: making y inout allows for already existing active
548     ! variables to become the target of a derivative assignment
549     !
550     subroutine setderiv_a0_a0(y,x)
551     type(active), intent(inout) :: y
552     type(active), intent(in) :: x
553     VECTOR_LOOP_VAR
554     VECTOR_LOOP_BEGIN
555     y%DELEM = x%DELEM
556     VECTOR_LOOP_END
557     end subroutine
558     subroutine setderiv_a1_a0(y,x)
559     type(active), intent(inout), dimension(:) :: y
560     type(active), intent(in) :: x
561     VECTOR_LOOP_VAR
562     VECTOR_LOOP_BEGIN
563     y%DELEM = x%DELEM
564     VECTOR_LOOP_END
565     end subroutine
566     subroutine setderiv_a1_a1(y,x)
567     type(active), intent(inout), dimension(:) :: y
568     type(active), intent(in), dimension(:) :: x
569     VECTOR_LOOP_VAR
570     VECTOR_LOOP_BEGIN
571     y%DELEM = x%DELEM
572     VECTOR_LOOP_END
573     end subroutine
574     subroutine setderiv_a2_a0(y,x)
575     type(active), intent(inout), dimension(:,:) :: y
576     type(active), intent(in) :: x
577     VECTOR_LOOP_VAR
578     VECTOR_LOOP_BEGIN
579     y%DELEM = x%DELEM
580     VECTOR_LOOP_END
581     end subroutine
582     subroutine setderiv_a2_a2(y,x)
583     type(active), intent(inout), dimension(:,:) :: y
584     type(active), intent(in), dimension(:,:) :: x
585     VECTOR_LOOP_VAR
586     VECTOR_LOOP_BEGIN
587     y%DELEM = x%DELEM
588     VECTOR_LOOP_END
589     end subroutine
590     subroutine setderiv_a3_a3(y,x)
591     type(active), intent(inout), dimension(:,:,:) :: y
592     type(active), intent(in), dimension(:,:,:) :: x
593     VECTOR_LOOP_VAR
594     VECTOR_LOOP_BEGIN
595     y%DELEM = x%DELEM
596     VECTOR_LOOP_END
597     end subroutine
598     !
599     ! set the derivative of y to be the negated derivative of x
600     ! note: making y inout allows for already existing active
601     ! variables to become the target of a derivative assignment
602     !
603    
604     subroutine set_neg_deriv_a0_a0(y,x)
605     type(active), intent(inout) :: y
606     type(active), intent(in) :: x
607     VECTOR_LOOP_VAR
608     VECTOR_LOOP_BEGIN
609     y%DELEM = -x%DELEM
610     VECTOR_LOOP_END
611     end subroutine
612    
613     subroutine set_neg_deriv_a1_a1(y,x)
614     type(active), intent(inout), dimension(:) :: y
615     type(active), intent(in), dimension(:) :: x
616     VECTOR_LOOP_VAR
617     VECTOR_LOOP_BEGIN
618     y%DELEM = -x%DELEM
619     VECTOR_LOOP_END
620     end subroutine
621    
622     !
623     ! increment the derivative of y by the derivative of x
624     ! note: making y inout allows for already existing active
625     ! variables to become the target of a derivative assignment
626     !
627    
628     subroutine inc_deriv_a0_a0(y,x)
629     type(active), intent(inout) :: y
630     type(active), intent(in) :: x
631     VECTOR_LOOP_VAR
632     VECTOR_LOOP_BEGIN
633     y%DELEM=y%DELEM + x%DELEM
634     VECTOR_LOOP_END
635     end subroutine
636    
637     subroutine inc_deriv_a1_a1(y,x)
638     type(active), intent(inout), dimension(:) :: y
639     type(active), intent(in), dimension(:) :: x
640     VECTOR_LOOP_VAR
641     VECTOR_LOOP_BEGIN
642     y%DELEM=y%DELEM + x%DELEM
643     VECTOR_LOOP_END
644     end subroutine
645    
646     subroutine inc_deriv_a2_a2(y,x)
647     type(active), intent(inout), dimension(:,:) :: y
648     type(active), intent(in), dimension(:,:) :: x
649     VECTOR_LOOP_VAR
650     VECTOR_LOOP_BEGIN
651     y%DELEM=y%DELEM + x%DELEM
652     VECTOR_LOOP_END
653     end subroutine
654    
655     !
656     ! decrement the derivative of y by the derivative of x
657     ! note: making y inout allows for already existing active
658     ! variables to become the target of a derivative assignment
659     !
660    
661     subroutine dec_deriv_a0_a0(y,x)
662     type(active), intent(inout) :: y
663     type(active), intent(in) :: x
664     VECTOR_LOOP_VAR
665     VECTOR_LOOP_BEGIN
666     y%DELEM=y%DELEM - x%DELEM
667     VECTOR_LOOP_END
668     end subroutine
669    
670     subroutine dec_deriv_a1_a1(y,x)
671     type(active), intent(inout), dimension(:) :: y
672     type(active), intent(in), dimension(:) :: x
673     VECTOR_LOOP_VAR
674     VECTOR_LOOP_BEGIN
675     y%DELEM=y%DELEM - x%DELEM
676     VECTOR_LOOP_END
677     end subroutine
678    
679     subroutine dec_deriv_a2_a2(y,x)
680     type(active), intent(inout), dimension(:,:) :: y
681     type(active), intent(in), dimension(:,:) :: x
682     VECTOR_LOOP_VAR
683     VECTOR_LOOP_BEGIN
684     y%DELEM=y%DELEM - x%DELEM
685     VECTOR_LOOP_END
686     end subroutine
687    
688     !
689     ! set derivative components to 0.0
690     !
691     subroutine zero_deriv_a0(x)
692     type(active), intent(inout) :: x
693     VECTOR_LOOP_VAR
694     VECTOR_LOOP_BEGIN
695     x%DELEM=0.0d0
696     VECTOR_LOOP_END
697     end subroutine
698    
699     subroutine zero_deriv_a1(x)
700     type(active), dimension(:), intent(inout) :: x
701     VECTOR_LOOP_VAR
702     VECTOR_LOOP_BEGIN
703     x%DELEM=0.0d0
704     VECTOR_LOOP_END
705     end subroutine
706    
707     subroutine zero_deriv_a2(x)
708     type(active), dimension(:,:), intent(inout) :: x
709     VECTOR_LOOP_VAR
710     VECTOR_LOOP_BEGIN
711     x%DELEM=0.0d0
712     VECTOR_LOOP_END
713     end subroutine
714    
715     subroutine zero_deriv_a3(x)
716     type(active), dimension(:,:,:), intent(inout) :: x
717     VECTOR_LOOP_VAR
718     VECTOR_LOOP_BEGIN
719     x%DELEM=0.0d0
720     VECTOR_LOOP_END
721     end subroutine
722    
723     subroutine zero_deriv_a4(x)
724     type(active), dimension(:,:,:,:), intent(inout) :: x
725     VECTOR_LOOP_VAR
726     VECTOR_LOOP_BEGIN
727     x%DELEM=0.0d0
728     VECTOR_LOOP_END
729     end subroutine
730    
731     #endif
732     !
733     ! conversions
734     !
735     subroutine convert_d0_a0(convertTo, convertFrom)
736     real(w2f__8), intent(out) :: convertTo
737     type(active), intent(in) :: convertFrom
738     convertTo=convertFrom%v
739     end subroutine
740     subroutine convert_d1_a1(convertTo, convertFrom)
741     real(w2f__8), dimension(:), intent(out) :: convertTo
742     type(active), dimension(:), intent(in) :: convertFrom
743     convertTo=convertFrom%v
744     end subroutine
745     subroutine convert_d2_a2(convertTo, convertFrom)
746     real(w2f__8), dimension(:,:), intent(out) :: convertTo
747     type(active), dimension(:,:), intent(in) :: convertFrom
748     convertTo=convertFrom%v
749     end subroutine
750     subroutine convert_d3_a3(convertTo, convertFrom)
751     real(w2f__8), dimension(:,:,:), intent(out) :: convertTo
752     type(active), dimension(:,:,:), intent(in) :: convertFrom
753     convertTo=convertFrom%v
754     end subroutine
755     subroutine convert_d4_a4(convertTo, convertFrom)
756     real(w2f__8), dimension(:,:,:,:), intent(out) :: convertTo
757     type(active), dimension(:,:,:,:), intent(in) :: convertFrom
758     convertTo=convertFrom%v
759     end subroutine
760     subroutine convert_d5_a5(convertTo, convertFrom)
761     real(w2f__8), dimension(:,:,:,:,:), intent(out) :: convertTo
762     type(active), dimension(:,:,:,:,:), intent(in) :: convertFrom
763     convertTo=convertFrom%v
764     end subroutine
765     subroutine convert_d6_a6(convertTo, convertFrom)
766     real(w2f__8), dimension(:,:,:,:,:,:), intent(out) :: convertTo
767     type(active), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
768     convertTo=convertFrom%v
769     end subroutine
770     subroutine convert_d7_a7(convertTo, convertFrom)
771     real(w2f__8), dimension(:,:,:,:,:,:,:), intent(out) :: convertTo
772     type(active), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
773     convertTo=convertFrom%v
774     end subroutine
775    
776     subroutine convert_a0_d0(convertTo, convertFrom)
777     type(active), intent(inout) :: convertTo
778     real(w2f__8), intent(in) :: convertFrom
779     convertTo%v=convertFrom
780     end subroutine
781     subroutine convert_a1_d1(convertTo, convertFrom)
782     type(active), dimension(:), intent(inout) :: convertTo
783     real(w2f__8), dimension(:), intent(in) :: convertFrom
784     convertTo%v=convertFrom
785     end subroutine
786     subroutine convert_a2_d2(convertTo, convertFrom)
787     type(active), dimension(:,:), intent(inout) :: convertTo
788     real(w2f__8), dimension(:,:), intent(in) :: convertFrom
789     convertTo%v=convertFrom
790     end subroutine
791     subroutine convert_a3_d3(convertTo, convertFrom)
792     type(active), dimension(:,:,:), intent(inout) :: convertTo
793     real(w2f__8), dimension(:,:,:), intent(in) :: convertFrom
794     convertTo%v=convertFrom
795     end subroutine
796     subroutine convert_a4_d4(convertTo, convertFrom)
797     type(active), dimension(:,:,:,:), intent(inout) :: convertTo
798     real(w2f__8), dimension(:,:,:,:), intent(in) :: convertFrom
799     convertTo%v=convertFrom
800     end subroutine
801     subroutine convert_a5_d5(convertTo, convertFrom)
802     type(active), dimension(:,:,:,:,:), intent(inout) :: convertTo
803     real(w2f__8), dimension(:,:,:,:,:), intent(in) :: convertFrom
804     convertTo%v=convertFrom
805     end subroutine
806     subroutine convert_a6_d6(convertTo, convertFrom)
807     type(active), dimension(:,:,:,:,:,:), intent(inout) :: convertTo
808     real(w2f__8), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
809     convertTo%v=convertFrom
810     end subroutine
811     subroutine convert_a7_d7(convertTo, convertFrom)
812     type(active), dimension(:,:,:,:,:,:,:), intent(inout) :: convertTo
813     real(w2f__8), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
814     convertTo%v=convertFrom
815     end subroutine
816     #ifndef DEFAULT_R8
817     subroutine convert_r0_a0(convertTo, convertFrom)
818     real(w2f__4), intent(out) :: convertTo
819     type(active), intent(in) :: convertFrom
820     convertTo=convertFrom%v
821     end subroutine
822     subroutine convert_r1_a1(convertTo, convertFrom)
823     real(w2f__4), dimension(:), intent(out) :: convertTo
824     type(active), dimension(:), intent(in) :: convertFrom
825     convertTo=convertFrom%v
826     end subroutine
827     subroutine convert_r2_a2(convertTo, convertFrom)
828     real(w2f__4), dimension(:,:), intent(out) :: convertTo
829     type(active), dimension(:,:), intent(in) :: convertFrom
830     convertTo=convertFrom%v
831     end subroutine
832     subroutine convert_r3_a3(convertTo, convertFrom)
833     real(w2f__4), dimension(:,:,:), intent(out) :: convertTo
834     type(active), dimension(:,:,:), intent(in) :: convertFrom
835     convertTo=convertFrom%v
836     end subroutine
837     subroutine convert_r4_a4(convertTo, convertFrom)
838     real(w2f__4), dimension(:,:,:,:), intent(out) :: convertTo
839     type(active), dimension(:,:,:,:), intent(in) :: convertFrom
840     convertTo=convertFrom%v
841     end subroutine
842     subroutine convert_r5_a5(convertTo, convertFrom)
843     real(w2f__4), dimension(:,:,:,:,:), intent(out) :: convertTo
844     type(active), dimension(:,:,:,:,:), intent(in) :: convertFrom
845     convertTo=convertFrom%v
846     end subroutine
847     subroutine convert_r6_a6(convertTo, convertFrom)
848     real(w2f__4), dimension(:,:,:,:,:,:), intent(out) :: convertTo
849     type(active), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
850     convertTo=convertFrom%v
851     end subroutine
852     subroutine convert_r7_a7(convertTo, convertFrom)
853     real(w2f__4), dimension(:,:,:,:,:,:,:), intent(out) :: convertTo
854     type(active), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
855     convertTo=convertFrom%v
856     end subroutine
857    
858     subroutine convert_a0_r0(convertTo, convertFrom)
859     type(active), intent(inout) :: convertTo
860     real(w2f__4), intent(in) :: convertFrom
861     convertTo%v=convertFrom
862     end subroutine
863     subroutine convert_a1_r1(convertTo, convertFrom)
864     type(active), dimension(:), intent(inout) :: convertTo
865     real(w2f__4), dimension(:), intent(in) :: convertFrom
866     convertTo%v=convertFrom
867     end subroutine
868     subroutine convert_a2_r2(convertTo, convertFrom)
869     type(active), dimension(:,:), intent(inout) :: convertTo
870     real(w2f__4), dimension(:,:), intent(in) :: convertFrom
871     convertTo%v=convertFrom
872     end subroutine
873     subroutine convert_a3_r3(convertTo, convertFrom)
874     type(active), dimension(:,:,:), intent(inout) :: convertTo
875     real(w2f__4), dimension(:,:,:), intent(in) :: convertFrom
876     convertTo%v=convertFrom
877     end subroutine
878     subroutine convert_a4_r4(convertTo, convertFrom)
879     type(active), dimension(:,:,:,:), intent(inout) :: convertTo
880     real(w2f__4), dimension(:,:,:,:), intent(in) :: convertFrom
881     convertTo%v=convertFrom
882     end subroutine
883     subroutine convert_a5_r5(convertTo, convertFrom)
884     type(active), dimension(:,:,:,:,:), intent(inout) :: convertTo
885     real(w2f__4), dimension(:,:,:,:,:), intent(in) :: convertFrom
886     convertTo%v=convertFrom
887     end subroutine
888     subroutine convert_a6_r6(convertTo, convertFrom)
889     type(active), dimension(:,:,:,:,:,:), intent(inout) :: convertTo
890     real(w2f__4), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
891     convertTo%v=convertFrom
892     end subroutine
893     subroutine convert_a7_r7(convertTo, convertFrom)
894     type(active), dimension(:,:,:,:,:,:,:), intent(inout) :: convertTo
895     real(w2f__4), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
896     convertTo%v=convertFrom
897     end subroutine
898     #endif
899     !
900     ! allocations
901     !
902     subroutine allocateMatching_d1_d1(toBeAllocated,allocateMatching)
903     implicit none
904     real(w2f__8), dimension(:), allocatable :: toBeAllocated
905     real(w2f__8), dimension(:) :: allocateMatching
906     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
907     allocate(toBeAllocated(size(allocateMatching)))
908     end subroutine
909     subroutine allocateMatching_a1_d1(toBeAllocated,allocateMatching)
910     implicit none
911     type(active), dimension(:), allocatable :: toBeAllocated
912     real(w2f__8), dimension(:) :: allocateMatching
913     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
914     allocate(toBeAllocated(size(allocateMatching)))
915     end subroutine
916     subroutine allocateMatching_d1_a1(toBeAllocated,allocateMatching)
917     implicit none
918     real(w2f__8), dimension(:), allocatable :: toBeAllocated
919     type(active), dimension(:) :: allocateMatching
920     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
921     allocate(toBeAllocated(size(allocateMatching)))
922     end subroutine
923     subroutine allocateMatching_a1_a1(toBeAllocated,allocateMatching)
924     implicit none
925     type(active), dimension(:), allocatable :: toBeAllocated
926     type(active), dimension(:) :: allocateMatching
927     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
928     allocate(toBeAllocated(size(allocateMatching)))
929     end subroutine
930     subroutine allocateMatching_d2_d2(toBeAllocated,allocateMatching)
931     implicit none
932     real(w2f__8), dimension(:,:), allocatable :: toBeAllocated
933     real(w2f__8), dimension(:,:) :: allocateMatching
934     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
935     allocate(toBeAllocated(size(allocateMatching,1), &
936     size(allocateMatching,2)))
937     end subroutine
938     subroutine allocateMatching_a2_d2(toBeAllocated,allocateMatching)
939     implicit none
940     type(active), dimension(:,:), allocatable :: toBeAllocated
941     real(w2f__8), dimension(:,:) :: allocateMatching
942     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
943     allocate(toBeAllocated(size(allocateMatching,1), &
944     size(allocateMatching,2)))
945     end subroutine
946     subroutine allocateMatching_d2_a2(toBeAllocated,allocateMatching)
947     implicit none
948     real(w2f__8), dimension(:,:), allocatable :: toBeAllocated
949     type(active), dimension(:,:) :: allocateMatching
950     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
951     allocate(toBeAllocated(size(allocateMatching,1), &
952     size(allocateMatching,2)))
953     end subroutine
954     subroutine allocateMatching_a2_a2(toBeAllocated,allocateMatching)
955     implicit none
956     type(active), dimension(:,:), allocatable :: toBeAllocated
957     type(active), dimension(:,:) :: allocateMatching
958     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
959     allocate(toBeAllocated(size(allocateMatching,1), &
960     size(allocateMatching,2)))
961     end subroutine
962     subroutine allocateMatching_d3_d3(toBeAllocated,allocateMatching)
963     implicit none
964     real(w2f__8), dimension(:,:,:), allocatable :: toBeAllocated
965     real(w2f__8), dimension(:,:,:) :: allocateMatching
966     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
967     allocate(toBeAllocated(size(allocateMatching,1), &
968     size(allocateMatching,2),&
969     size(allocateMatching,3)))
970     end subroutine
971     subroutine allocateMatching_a3_d3(toBeAllocated,allocateMatching)
972     implicit none
973     type(active), dimension(:,:,:), allocatable :: toBeAllocated
974     real(w2f__8), dimension(:,:,:) :: allocateMatching
975     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
976     allocate(toBeAllocated(size(allocateMatching,1), &
977     size(allocateMatching,2),&
978     size(allocateMatching,3)))
979     end subroutine
980     subroutine allocateMatching_d3_a3(toBeAllocated,allocateMatching)
981     implicit none
982     real(w2f__8), dimension(:,:,:), allocatable :: toBeAllocated
983     type(active), dimension(:,:,:) :: allocateMatching
984     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
985     allocate(toBeAllocated(size(allocateMatching,1), &
986     size(allocateMatching,2),&
987     size(allocateMatching,3)))
988     end subroutine
989     subroutine allocateMatching_a3_a3(toBeAllocated,allocateMatching)
990     implicit none
991     type(active), dimension(:,:,:), allocatable :: toBeAllocated
992     type(active), dimension(:,:,:) :: allocateMatching
993     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
994     allocate(toBeAllocated(size(allocateMatching,1), &
995     size(allocateMatching,2),&
996     size(allocateMatching,3)))
997     end subroutine
998     subroutine allocateMatching_a4_a4(toBeAllocated,allocateMatching)
999     implicit none
1000     type(active), dimension(:,:,:,:), allocatable :: toBeAllocated
1001     type(active), dimension(:,:,:,:) :: allocateMatching
1002     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1003     allocate(toBeAllocated(size(allocateMatching,1), &
1004     size(allocateMatching,2),&
1005     size(allocateMatching,3),&
1006     size(allocateMatching,4)))
1007     end subroutine
1008     subroutine allocateMatching_d4_a4(toBeAllocated,allocateMatching)
1009     implicit none
1010     real(w2f__8), dimension(:,:,:,:), allocatable :: toBeAllocated
1011     type(active), dimension(:,:,:,:) :: allocateMatching
1012     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1013     allocate(toBeAllocated(size(allocateMatching,1), &
1014     size(allocateMatching,2),&
1015     size(allocateMatching,3),&
1016     size(allocateMatching,4)))
1017     end subroutine
1018     subroutine allocateMatching_d5_d5(toBeAllocated,allocateMatching)
1019     implicit none
1020     real(w2f__8), dimension(:,:,:,:,:), allocatable :: toBeAllocated
1021     real(w2f__8), dimension(:,:,:,:,:) :: allocateMatching
1022     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1023     allocate(toBeAllocated(size(allocateMatching,1), &
1024     size(allocateMatching,2),&
1025     size(allocateMatching,3),&
1026     size(allocateMatching,4),&
1027     size(allocateMatching,5)))
1028     end subroutine
1029     subroutine allocateMatching_a5_d5(toBeAllocated,allocateMatching)
1030     implicit none
1031     type(active), dimension(:,:,:,:,:), allocatable :: toBeAllocated
1032     real(w2f__8), dimension(:,:,:,:,:) :: allocateMatching
1033     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1034     allocate(toBeAllocated(size(allocateMatching,1), &
1035     size(allocateMatching,2),&
1036     size(allocateMatching,3),&
1037     size(allocateMatching,4),&
1038     size(allocateMatching,5)))
1039     end subroutine
1040     subroutine allocateMatching_d5_a5(toBeAllocated,allocateMatching)
1041     implicit none
1042     real(w2f__8), dimension(:,:,:,:,:), allocatable :: toBeAllocated
1043     type(active), dimension(:,:,:,:,:) :: allocateMatching
1044     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1045     allocate(toBeAllocated(size(allocateMatching,1), &
1046     size(allocateMatching,2),&
1047     size(allocateMatching,3),&
1048     size(allocateMatching,4),&
1049     size(allocateMatching,5)))
1050     end subroutine
1051     subroutine allocateMatching_a5_a5(toBeAllocated,allocateMatching)
1052     implicit none
1053     type(active), dimension(:,:,:,:,:), allocatable :: toBeAllocated
1054     type(active), dimension(:,:,:,:,:) :: allocateMatching
1055     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1056     allocate(toBeAllocated(size(allocateMatching,1), &
1057     size(allocateMatching,2),&
1058     size(allocateMatching,3),&
1059     size(allocateMatching,4),&
1060     size(allocateMatching,5)))
1061     end subroutine
1062     subroutine allocateMatching_d6_d6(toBeAllocated,allocateMatching)
1063     implicit none
1064     real(w2f__8), dimension(:,:,:,:,:,:), allocatable :: toBeAllocated
1065     real(w2f__8), dimension(:,:,:,:,:,:) :: allocateMatching
1066     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1067     allocate(toBeAllocated(size(allocateMatching,1), &
1068     size(allocateMatching,2),&
1069     size(allocateMatching,3),&
1070     size(allocateMatching,4),&
1071     size(allocateMatching,5),&
1072     size(allocateMatching,6)))
1073     end subroutine
1074     subroutine allocateMatching_a6_d6(toBeAllocated,allocateMatching)
1075     implicit none
1076     type(active), dimension(:,:,:,:,:,:), allocatable :: toBeAllocated
1077     real(w2f__8), dimension(:,:,:,:,:,:) :: allocateMatching
1078     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1079     allocate(toBeAllocated(size(allocateMatching,1), &
1080     size(allocateMatching,2),&
1081     size(allocateMatching,3),&
1082     size(allocateMatching,4),&
1083     size(allocateMatching,5),&
1084     size(allocateMatching,6)))
1085     end subroutine
1086     subroutine allocateMatching_d6_a6(toBeAllocated,allocateMatching)
1087     implicit none
1088     real(w2f__8), dimension(:,:,:,:,:,:), allocatable :: toBeAllocated
1089     type(active), dimension(:,:,:,:,:,:) :: allocateMatching
1090     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1091     allocate(toBeAllocated(size(allocateMatching,1), &
1092     size(allocateMatching,2),&
1093     size(allocateMatching,3),&
1094     size(allocateMatching,4),&
1095     size(allocateMatching,5),&
1096     size(allocateMatching,6)))
1097     end subroutine
1098     subroutine allocateMatching_a6_a6(toBeAllocated,allocateMatching)
1099     implicit none
1100     type(active), dimension(:,:,:,:,:,:), allocatable :: toBeAllocated
1101     type(active), dimension(:,:,:,:,:,:) :: allocateMatching
1102     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1103     allocate(toBeAllocated(size(allocateMatching,1), &
1104     size(allocateMatching,2),&
1105     size(allocateMatching,3),&
1106     size(allocateMatching,4),&
1107     size(allocateMatching,5),&
1108     size(allocateMatching,6)))
1109     end subroutine
1110     #ifndef DEFAULT_R8
1111     subroutine allocateMatching_r1_r1(toBeAllocated,allocateMatching)
1112     implicit none
1113     real(w2f__4), dimension(:), allocatable :: toBeAllocated
1114     real(w2f__4), dimension(:) :: allocateMatching
1115     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1116     allocate(toBeAllocated(size(allocateMatching)))
1117     end subroutine
1118     subroutine allocateMatching_d1_r1(toBeAllocated,allocateMatching)
1119     implicit none
1120     real(w2f__8), dimension(:), allocatable :: toBeAllocated
1121     real(w2f__4), dimension(:) :: allocateMatching
1122     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1123     allocate(toBeAllocated(size(allocateMatching)))
1124     end subroutine
1125     subroutine allocateMatching_r1_d1(toBeAllocated,allocateMatching)
1126     implicit none
1127     real(w2f__4), dimension(:), allocatable :: toBeAllocated
1128     real(w2f__8), dimension(:) :: allocateMatching
1129     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1130     allocate(toBeAllocated(size(allocateMatching)))
1131     end subroutine
1132     subroutine allocateMatching_a1_r1(toBeAllocated,allocateMatching)
1133     implicit none
1134     type(active), dimension(:), allocatable :: toBeAllocated
1135     real(w2f__4), dimension(:) :: allocateMatching
1136     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1137     allocate(toBeAllocated(size(allocateMatching)))
1138     end subroutine
1139     subroutine allocateMatching_r1_a1(toBeAllocated,allocateMatching)
1140     implicit none
1141     real(w2f__4), dimension(:), allocatable :: toBeAllocated
1142     type(active), dimension(:) :: allocateMatching
1143     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1144     allocate(toBeAllocated(size(allocateMatching)))
1145     end subroutine
1146     subroutine allocateMatching_a2_r2(toBeAllocated,allocateMatching)
1147     implicit none
1148     type(active), dimension(:,:), allocatable :: toBeAllocated
1149     real(w2f__4), dimension(:,:) :: allocateMatching
1150     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1151     allocate(toBeAllocated(size(allocateMatching,1), &
1152     size(allocateMatching,2)))
1153     end subroutine
1154     subroutine allocateMatching_r2_a2(toBeAllocated,allocateMatching)
1155     implicit none
1156     real(w2f__4), dimension(:,:), allocatable :: toBeAllocated
1157     type(active), dimension(:,:) :: allocateMatching
1158     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1159     allocate(toBeAllocated(size(allocateMatching,1), &
1160     size(allocateMatching,2)))
1161     end subroutine
1162     #endif
1163     !
1164     ! allocate shape
1165     !
1166     subroutine allocateShape_d1(toBeAllocated,s1)
1167     implicit none
1168     real(w2f__8), dimension(:), allocatable :: toBeAllocated
1169     integer(w2f__i8) :: s1
1170     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1171     allocate(toBeAllocated(s1))
1172     end subroutine
1173     subroutine allocateShape_d2(toBeAllocated,s1,s2)
1174     implicit none
1175     real(w2f__8), dimension(:,:), allocatable :: toBeAllocated
1176     integer(w2f__i8) :: s1,s2
1177     if (allocated(toBeAllocated)) deallocate(toBeAllocated)
1178     allocate(toBeAllocated(s1,s2))
1179     end subroutine
1180     !
1181     ! shape tests
1182     !
1183     subroutine shapeTest_a1_d1(allocatedVar,origVar)
1184     implicit none
1185     type(active), dimension(:), allocatable :: allocatedVar
1186     real(w2f__8), dimension(:) :: origVar
1187     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1188     end subroutine
1189     subroutine shapeTest_a1_a1(allocatedVar,origVar)
1190     implicit none
1191     type(active), dimension(:), allocatable :: allocatedVar
1192     type(active), dimension(:) :: origVar
1193     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1194     end subroutine
1195     subroutine shapeTest_d1_a1(allocatedVar,origVar)
1196     implicit none
1197     real(w2f__8), dimension(:), allocatable :: allocatedVar
1198     type(active), dimension(:) :: origVar
1199     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1200     end subroutine
1201     subroutine shapeTest_a2_d2(allocatedVar,origVar)
1202     implicit none
1203     type(active), dimension(:,:), allocatable :: allocatedVar
1204     real(w2f__8), dimension(:,:) :: origVar
1205     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1206     end subroutine
1207     subroutine shapeTest_a2_a2(allocatedVar,origVar)
1208     implicit none
1209     type(active), dimension(:,:), allocatable :: allocatedVar
1210     type(active), dimension(:,:) :: origVar
1211     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1212     end subroutine
1213     subroutine shapeTest_d2_a2(allocatedVar,origVar)
1214     implicit none
1215     real(w2f__8), dimension(:,:), allocatable :: allocatedVar
1216     type(active), dimension(:,:) :: origVar
1217     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1218     end subroutine
1219     subroutine shapeTest_a3_a3(allocatedVar,origVar)
1220     implicit none
1221     type(active), dimension(:,:,:), allocatable :: allocatedVar
1222     type(active), dimension(:,:,:) :: origVar
1223     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1224     end subroutine
1225     subroutine shapeTest_d3_a3(allocatedVar,origVar)
1226     implicit none
1227     real(w2f__8), dimension(:,:,:), allocatable :: allocatedVar
1228     type(active), dimension(:,:,:) :: origVar
1229     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1230     end subroutine
1231     subroutine shapeTest_a4_a4(allocatedVar,origVar)
1232     implicit none
1233     type(active), dimension(:,:,:,:), allocatable :: allocatedVar
1234     type(active), dimension(:,:,:,:) :: origVar
1235     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1236     end subroutine
1237     subroutine shapeTest_d4_a4(allocatedVar,origVar)
1238     implicit none
1239     real(w2f__8), dimension(:,:,:,:), allocatable :: allocatedVar
1240     type(active), dimension(:,:,:,:) :: origVar
1241     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1242     end subroutine
1243     subroutine shapeTest_a5_a5(allocatedVar,origVar)
1244     implicit none
1245     type(active), dimension(:,:,:,:,:), allocatable :: allocatedVar
1246     type(active), dimension(:,:,:,:,:) :: origVar
1247     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1248     end subroutine
1249     subroutine shapeTest_d5_a5(allocatedVar,origVar)
1250     implicit none
1251     real(w2f__8), dimension(:,:,:,:,:), allocatable :: allocatedVar
1252     type(active), dimension(:,:,:,:,:) :: origVar
1253     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1254     end subroutine
1255     subroutine shapeTest_a5_d5(allocatedVar,origVar)
1256     implicit none
1257     type(active), dimension(:,:,:,:,:), allocatable :: allocatedVar
1258     real(w2f__8), dimension(:,:,:,:,:) :: origVar
1259     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1260     end subroutine
1261     subroutine shapeTest_a6_a6(allocatedVar,origVar)
1262     implicit none
1263     type(active), dimension(:,:,:,:,:,:), allocatable :: allocatedVar
1264     type(active), dimension(:,:,:,:,:,:) :: origVar
1265     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1266     end subroutine
1267     subroutine shapeTest_d6_a6(allocatedVar,origVar)
1268     implicit none
1269     real(w2f__8), dimension(:,:,:,:,:,:), allocatable :: allocatedVar
1270     type(active), dimension(:,:,:,:,:,:) :: origVar
1271     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1272     end subroutine
1273     subroutine shapeTest_a6_d6(allocatedVar,origVar)
1274     implicit none
1275     type(active), dimension(:,:,:,:,:,:), allocatable :: allocatedVar
1276     real(w2f__8), dimension(:,:,:,:,:,:) :: origVar
1277     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1278     end subroutine
1279     #ifndef DEFAULT_R8
1280     subroutine shapeTest_r1_a1(allocatedVar,origVar)
1281     implicit none
1282     real(w2f__4), dimension(:), allocatable :: allocatedVar
1283     type(active), dimension(:) :: origVar
1284     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1285     end subroutine
1286     subroutine shapeTest_r2_a2(allocatedVar,origVar)
1287     implicit none
1288     real(w2f__4), dimension(:,:), allocatable :: allocatedVar
1289     type(active), dimension(:,:) :: origVar
1290     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1291     end subroutine
1292     subroutine shapeTest_a1_r1(allocatedVar,origVar)
1293     implicit none
1294     type(active), dimension(:), allocatable :: allocatedVar
1295     real(w2f__4), dimension(:) :: origVar
1296     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1297     end subroutine
1298     subroutine shapeTest_a2_r2(allocatedVar,origVar)
1299     implicit none
1300     type(active), dimension(:,:), allocatable :: allocatedVar
1301     real(w2f__4), dimension(:,:) :: origVar
1302     if (.not. all(shape(allocatedVar)==shape(origVar))) call runTimeErrorStop(shapeChange)
1303     end subroutine
1304     #endif
1305     subroutine runTimeErrorStopI(mesgId)
1306     implicit none
1307     integer mesgId
1308     select case (mesgId)
1309     case (shapeChange)
1310     stop "ERROR: OAD run time library: detected shape change"
1311     end select
1312     end subroutine
1313    
1314     end module

  ViewVC Help
Powered by ViewVC 1.1.22