Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:45:09 UTC

view on githubraw file Latest commit 8702af1f on 2012-09-20 23:12:48 UTC
8702af1f36 Patr*0001 !#########################################################
                0002 ! This file is part of OpenAD released under the LGPL.   #
                0003 ! The full COPYRIGHT notice can be found in the top      #
                0004 ! level directory of the OpenAD distribution             #
                0005 !#########################################################
                0006         module OAD_active
                0007         use w2f__types
                0008         implicit none
                0009         private :: runTimeErrorStop, shapeChange
                0010         public :: active
                0011 #ifndef TRACE
                0012         public :: saxpy, sax, zero_deriv, setderiv
                0013         public :: set_neg_deriv, inc_deriv, dec_deriv
                0014 #endif
                0015         public :: oad_convert
                0016         public :: oad_allocateMatching, oad_allocateShape, oad_shapeTest
                0017 #ifndef TRACE
                0018         integer :: count_mult = 0
                0019         integer :: count_add = 0
                0020 #endif
                0021         integer, parameter :: shapeChange=0
                0022 #ifdef VECTOR
                0023         integer :: max_deriv_vec_len
                0024         parameter ( max_deriv_vec_len = 100 )
                0025         integer :: oad_ctmp_i1, oad_ctmp_i2, oad_ctmp_i3, & 
                0026         oad_ctmp_i4, oad_ctmp_i5, oad_ctmp_i6, oad_ctmp_i7
                0027 # define VECTOR_DIM , dimension(max_deriv_vec_len)
                0028 # define VECTOR_LOOP_VAR integer :: i
                0029 # define VECTOR_LOOP_BEGIN do i=1,max_deriv_vec_len
                0030 # define VECTOR_LOOP_END end do
                0031 # define DELEM d(i)
                0032 #else 
                0033 # define VECTOR_DIM 
                0034 # define VECTOR_LOOP_VAR
                0035 # define VECTOR_LOOP_BEGIN
                0036 # define VECTOR_LOOP_END
                0037 # define DELEM d
                0038 #endif
                0039 #ifdef SCALARNDI
                0040 # define DINIT
                0041 #else
                0042 # define DINIT =0.0d0
                0043 #endif
                0044 
                0045         !
                0046         ! active needs to be a sequence type
                0047         !  with no initialization
                0048         !
                0049         type active
                0050           sequence
                0051           real(w2f__8) :: v 
                0052 #ifndef TRACE
                0053           ! initialization does not work for active variables
                0054           ! inside of common block, such as in boxmodel
                0055           ! initialization is required for correct adjoint
                0056           real(w2f__8) VECTOR_DIM :: d DINIT
                0057 #endif
                0058         end type
                0059 #ifndef TRACE
                0060         interface saxpy
                0061           module procedure saxpy_d0_a0_a0
                0062           module procedure saxpy_l0_a0_a0
                0063           module procedure saxpy_i0_a0_a0
                0064           module procedure saxpy_d0_a0_a1
                0065           module procedure saxpy_d0_a1_a1
                0066           module procedure saxpy_d0_a2_a2
                0067           module procedure saxpy_d1_a0_a1
                0068           module procedure saxpy_d1_a1_a1 
                0069           module procedure saxpy_l1_a1_a1 
                0070           module procedure saxpy_i1_a1_a1
                0071           module procedure saxpy_a1_a1_a1
                0072           module procedure saxpy_d2_a0_a2
                0073           module procedure saxpy_d2_a2_a2
                0074 # ifndef DEFAULT_R8
                0075           module procedure saxpy_r0_a0_a0
                0076           module procedure saxpy_r0_a1_a1
                0077           module procedure saxpy_r1_a0_a1
                0078           module procedure saxpy_r1_a1_a1
                0079 # endif
                0080         end interface
                0081         
                0082         interface setderiv
                0083           module procedure setderiv_a0_a0
                0084           module procedure setderiv_a1_a0
                0085           module procedure setderiv_a1_a1
                0086           module procedure setderiv_a2_a0
                0087           module procedure setderiv_a2_a2
                0088           module procedure setderiv_a3_a3
                0089         end interface
                0090 
                0091         interface set_neg_deriv
                0092           module procedure set_neg_deriv_a0_a0
                0093           module procedure set_neg_deriv_a1_a1
                0094         end interface set_neg_deriv
                0095 
                0096         interface inc_deriv
                0097           module procedure inc_deriv_a0_a0
                0098           module procedure inc_deriv_a1_a1
                0099           module procedure inc_deriv_a2_a2
                0100         end interface inc_deriv
                0101 
                0102         interface dec_deriv
                0103           module procedure dec_deriv_a0_a0
                0104           module procedure dec_deriv_a1_a1
                0105           module procedure dec_deriv_a2_a2
                0106         end interface dec_deriv
                0107 
                0108         interface zero_deriv
                0109           module procedure zero_deriv_a0
                0110           module procedure zero_deriv_a1
                0111           module procedure zero_deriv_a2
                0112           module procedure zero_deriv_a3
                0113           module procedure zero_deriv_a4
                0114         end interface
                0115 
                0116         interface sax
                0117           module procedure sax_d0_a0_a0 
                0118           module procedure sax_l0_a0_a0
                0119           module procedure sax_i0_a0_a0
                0120           module procedure sax_d0_a0_a1
                0121           module procedure sax_d0_a0_a2
                0122           module procedure sax_d0_a1_a1
                0123           module procedure sax_d0_a2_a2
                0124           module procedure sax_d0_a3_a3
                0125           module procedure sax_d1_a0_a1
                0126           module procedure sax_d1_a1_a1 
                0127           module procedure sax_l1_a1_a1 
                0128           module procedure sax_i1_a1_a1
                0129           module procedure sax_d2_a0_a2
                0130           module procedure sax_d2_a2_a2
                0131 # ifndef DEFAULT_R8
                0132           module procedure sax_r0_a0_a0
                0133 # endif
                0134         end interface
                0135 
                0136 #endif
                0137         interface oad_convert
                0138           module procedure convert_d0_a0
                0139           module procedure convert_d1_a1
                0140           module procedure convert_d2_a2
                0141           module procedure convert_d3_a3
                0142           module procedure convert_d4_a4
                0143           module procedure convert_d5_a5
                0144           module procedure convert_d6_a6
                0145           module procedure convert_d7_a7
                0146           module procedure convert_a0_d0
                0147           module procedure convert_a1_d1
                0148           module procedure convert_a2_d2
                0149           module procedure convert_a3_d3
                0150           module procedure convert_a4_d4
                0151           module procedure convert_a5_d5
                0152           module procedure convert_a6_d6
                0153           module procedure convert_a7_d7
                0154 #ifndef DEFAULT_R8
                0155           module procedure convert_r0_a0
                0156           module procedure convert_r1_a1
                0157           module procedure convert_r2_a2
                0158           module procedure convert_r3_a3
                0159           module procedure convert_r4_a4
                0160           module procedure convert_r5_a5
                0161           module procedure convert_r6_a6
                0162           module procedure convert_r7_a7
                0163           module procedure convert_a0_r0
                0164           module procedure convert_a1_r1
                0165           module procedure convert_a2_r2
                0166           module procedure convert_a3_r3
                0167           module procedure convert_a4_r4
                0168           module procedure convert_a5_r5
                0169           module procedure convert_a6_r6
                0170           module procedure convert_a7_r7
                0171 #endif
                0172         end interface
                0173 
                0174         interface oad_allocateMatching
                0175           module procedure allocateMatching_d1_d1
                0176           module procedure allocateMatching_a1_d1
                0177           module procedure allocateMatching_d1_a1
                0178           module procedure allocateMatching_a1_a1
                0179           module procedure allocateMatching_d2_d2
                0180           module procedure allocateMatching_a2_d2
                0181           module procedure allocateMatching_d2_a2
                0182           module procedure allocateMatching_a2_a2
                0183           module procedure allocateMatching_d3_d3
                0184           module procedure allocateMatching_a3_d3
                0185           module procedure allocateMatching_d3_a3
                0186           module procedure allocateMatching_a3_a3
                0187           module procedure allocateMatching_a4_a4
                0188           module procedure allocateMatching_d4_a4
                0189           module procedure allocateMatching_d5_d5
                0190           module procedure allocateMatching_a5_d5
                0191           module procedure allocateMatching_d5_a5
                0192           module procedure allocateMatching_a5_a5
                0193           module procedure allocateMatching_d6_d6
                0194           module procedure allocateMatching_a6_d6
                0195           module procedure allocateMatching_d6_a6
                0196           module procedure allocateMatching_a6_a6
                0197 #ifndef DEFAULT_R8
                0198           module procedure allocateMatching_r1_r1
                0199           module procedure allocateMatching_d1_r1
                0200           module procedure allocateMatching_r1_d1
                0201           module procedure allocateMatching_a1_r1
                0202           module procedure allocateMatching_r1_a1
                0203           module procedure allocateMatching_a2_r2
                0204           module procedure allocateMatching_r2_a2
                0205 #endif
                0206         end interface 
                0207 
                0208         interface oad_allocateShape
                0209           module procedure allocateShape_d1
                0210           module procedure allocateShape_d2
                0211         end interface 
                0212 
                0213         interface oad_shapeTest
                0214           module procedure shapeTest_a1_d1
                0215           module procedure shapeTest_a1_a1
                0216           module procedure shapeTest_d1_a1
                0217           module procedure shapeTest_a2_d2
                0218           module procedure shapeTest_a2_a2
                0219           module procedure shapeTest_d2_a2
                0220           module procedure shapeTest_a3_a3
                0221           module procedure shapeTest_d3_a3
                0222           module procedure shapeTest_a4_a4
                0223           module procedure shapeTest_d4_a4
                0224           module procedure shapeTest_a5_a5
                0225           module procedure shapeTest_d5_a5
                0226           module procedure shapeTest_a5_d5
                0227           module procedure shapeTest_a6_a6
                0228           module procedure shapeTest_d6_a6
                0229           module procedure shapeTest_a6_d6
                0230 #ifndef DEFAULT_R8
                0231           module procedure shapeTest_r1_a1
                0232           module procedure shapeTest_r2_a2
                0233           module procedure shapeTest_a1_r1
                0234           module procedure shapeTest_a2_r2
                0235 #endif
                0236         end interface 
                0237 
                0238         interface runTimeErrorStop
                0239           module procedure runTimeErrorStopI
                0240         end interface 
                0241 
                0242         contains
                0243 #ifndef TRACE
                0244         !
                0245         ! chain rule saxpy to be used in forward and reverse modes
                0246         !
                0247         subroutine saxpy_d0_a0_a0(a,x,y)
                0248           real(w2f__8), intent(in) :: a
                0249           type(active), intent(in) :: x
                0250           type(active), intent(inout) :: y
                0251           VECTOR_LOOP_VAR
                0252           VECTOR_LOOP_BEGIN
                0253             y%DELEM=y%DELEM + x%DELEM*a
                0254           VECTOR_LOOP_END
                0255         end subroutine
                0256         subroutine saxpy_l0_a0_a0(a,x,y)
                0257           integer(kind=w2f__i8), intent(in) :: a
                0258           type(active), intent(in) :: x
                0259           type(active), intent(inout) :: y
                0260           VECTOR_LOOP_VAR
                0261           VECTOR_LOOP_BEGIN
                0262             y%DELEM=y%DELEM + x%DELEM*a
                0263           VECTOR_LOOP_END
                0264         end subroutine
                0265         subroutine saxpy_i0_a0_a0(a,x,y)
                0266           integer(kind=w2f__i4), intent(in) :: a
                0267           type(active), intent(in) :: x
                0268           type(active), intent(inout) :: y
                0269           VECTOR_LOOP_VAR
                0270           VECTOR_LOOP_BEGIN
                0271             y%DELEM=y%DELEM + x%DELEM*a
                0272           VECTOR_LOOP_END
                0273         end subroutine
                0274         subroutine saxpy_d0_a0_a1(a,x,y)
                0275           real(w2f__8), intent(in) :: a
                0276           type(active), intent(in) :: x
                0277           type(active), dimension(:), intent(inout) :: y
                0278           VECTOR_LOOP_VAR
                0279           VECTOR_LOOP_BEGIN
                0280             y%DELEM=y%DELEM + x%DELEM*a
                0281           VECTOR_LOOP_END
                0282         end subroutine
                0283         subroutine saxpy_d0_a1_a1(a,x,y)
                0284           real(w2f__8), intent(in) :: a
                0285           type(active), dimension(:), intent(in) :: x
                0286           type(active), dimension(:), intent(inout) :: y
                0287           VECTOR_LOOP_VAR
                0288           VECTOR_LOOP_BEGIN
                0289             y%DELEM=y%DELEM + x%DELEM*a
                0290           VECTOR_LOOP_END
                0291         end subroutine
                0292         subroutine saxpy_d0_a2_a2(a,x,y)
                0293           real(w2f__8), intent(in) :: a
                0294           type(active), dimension(:,:), intent(in) :: x
                0295           type(active), dimension(:,:), intent(inout) :: y
                0296           VECTOR_LOOP_VAR
                0297           VECTOR_LOOP_BEGIN
                0298             y%DELEM=y%DELEM + x%DELEM*a
                0299           VECTOR_LOOP_END
                0300         end subroutine
                0301         subroutine saxpy_d1_a0_a1(a,x,y)
                0302           real(w2f__8), dimension(:), intent(in) :: a
                0303           type(active), intent(in) :: x
                0304           type(active), dimension(:), intent(inout) :: y
                0305           VECTOR_LOOP_VAR
                0306           VECTOR_LOOP_BEGIN
                0307             y%DELEM=y%DELEM+x%DELEM*a
                0308           VECTOR_LOOP_END
                0309         end subroutine
                0310         subroutine saxpy_d1_a1_a1(a,x,y)
                0311           real(w2f__8), dimension(:), intent(in) :: a
                0312           type(active), dimension(:), intent(in) :: x
                0313           type(active), dimension(:), intent(inout) :: y
                0314           VECTOR_LOOP_VAR
                0315           VECTOR_LOOP_BEGIN
                0316             y%DELEM=y%DELEM+x%DELEM*a
                0317           VECTOR_LOOP_END
                0318         end subroutine
                0319         subroutine saxpy_l1_a1_a1(a,x,y)
                0320           integer(kind=w2f__i8), dimension(:), intent(in) :: a
                0321           type(active), dimension(:), intent(in) :: x
                0322           type(active), dimension(:), intent(inout) :: y
                0323           VECTOR_LOOP_VAR
                0324           VECTOR_LOOP_BEGIN
                0325             y%DELEM=y%DELEM+x%DELEM*a
                0326           VECTOR_LOOP_END
                0327         end subroutine
                0328         subroutine saxpy_i1_a1_a1(a,x,y)
                0329           integer(kind=w2f__i4), dimension(:), intent(in) :: a
                0330           type(active), dimension(:), intent(in) :: x
                0331           type(active), dimension(:), intent(inout) :: y
                0332           VECTOR_LOOP_VAR
                0333           VECTOR_LOOP_BEGIN
                0334             y%DELEM=y%DELEM+x%DELEM*a
                0335           VECTOR_LOOP_END
                0336         end subroutine
                0337         subroutine saxpy_a1_a1_a1(a,x,y)
                0338           type(active), dimension(:), intent(in) :: a
                0339           type(active), dimension(:), intent(in) :: x
                0340           type(active), dimension(:), intent(inout) :: y
                0341           VECTOR_LOOP_VAR
                0342           VECTOR_LOOP_BEGIN
                0343             y%DELEM=y%DELEM+x%DELEM*a%v
                0344           VECTOR_LOOP_END
                0345         end subroutine
                0346         subroutine saxpy_d2_a0_a2(a,x,y)
                0347           real(w2f__8), dimension(:,:), intent(in) :: a
                0348           type(active), intent(in) :: x
                0349           type(active), dimension(:,:), intent(inout) :: y
                0350           VECTOR_LOOP_VAR
                0351           VECTOR_LOOP_BEGIN
                0352             y%DELEM=y%DELEM + x%DELEM*a
                0353           VECTOR_LOOP_END
                0354         end subroutine
                0355         subroutine saxpy_d2_a2_a2(a,x,y)
                0356           real(w2f__8), dimension(:,:), intent(in) :: a
                0357           type(active), dimension(:,:), intent(in) :: x
                0358           type(active), dimension(:,:), intent(inout) :: y
                0359           VECTOR_LOOP_VAR
                0360           VECTOR_LOOP_BEGIN
                0361             y%DELEM=y%DELEM + x%DELEM*a
                0362           VECTOR_LOOP_END
                0363         end subroutine
                0364 # ifndef DEFAULT_R8
                0365         subroutine saxpy_r0_a0_a0(a,x,y)
                0366           real(w2f__4), intent(in) :: a
                0367           type(active), intent(in) :: x
                0368           type(active), intent(inout) :: y
                0369           VECTOR_LOOP_VAR
                0370           VECTOR_LOOP_BEGIN
                0371             y%DELEM=y%DELEM + x%DELEM*a
                0372           VECTOR_LOOP_END
                0373         end subroutine
                0374         subroutine saxpy_r0_a1_a1(a,x,y)
                0375           real(w2f__4), intent(in) :: a
                0376           type(active), dimension(:), intent(in) :: x
                0377           type(active), dimension(:), intent(inout) :: y
                0378           VECTOR_LOOP_VAR
                0379           VECTOR_LOOP_BEGIN
                0380             y%DELEM=y%DELEM + x%DELEM*a
                0381           VECTOR_LOOP_END
                0382         end subroutine
                0383         subroutine saxpy_r1_a0_a1(a,x,y)
                0384           real(w2f__4), dimension(:), intent(in) :: a
                0385           type(active), intent(in) :: x
                0386           type(active), dimension(:), intent(inout) :: y
                0387           VECTOR_LOOP_VAR
                0388           VECTOR_LOOP_BEGIN
                0389             y%DELEM=y%DELEM + x%DELEM*a
                0390           VECTOR_LOOP_END
                0391         end subroutine
                0392         subroutine saxpy_r1_a1_a1(a,x,y)
                0393           real(w2f__4), dimension(:), intent(in) :: a
                0394           type(active), dimension(:), intent(in) :: x
                0395           type(active), dimension(:), intent(inout) :: y
                0396           VECTOR_LOOP_VAR
                0397           VECTOR_LOOP_BEGIN
                0398             y%DELEM=y%DELEM+x%DELEM*a
                0399           VECTOR_LOOP_END
                0400         end subroutine
                0401 # endif
                0402         !
                0403         ! chain rule saxpy to be used in forward and reverse modes
                0404         ! derivative component of y is equal to zero initially
                0405         ! note: y needs to be inout as otherwise value component gets
                0406         ! zeroed out
                0407         !
                0408         subroutine sax_d0_a0_a0(a,x,y)
                0409           real(w2f__8), intent(in) :: a
                0410           type(active), intent(in) :: x
                0411           type(active), intent(inout) :: y
                0412           VECTOR_LOOP_VAR
                0413           VECTOR_LOOP_BEGIN
                0414             y%DELEM=x%DELEM*a
                0415           VECTOR_LOOP_END
                0416         end subroutine
                0417         subroutine sax_l0_a0_a0(a,x,y)
                0418           integer(kind=w2f__i8), intent(in) :: a
                0419           type(active), intent(in) :: x
                0420           type(active), intent(inout) :: y
                0421           VECTOR_LOOP_VAR
                0422           VECTOR_LOOP_BEGIN
                0423             y%DELEM=x%DELEM*a
                0424           VECTOR_LOOP_END
                0425         end subroutine
                0426         subroutine sax_i0_a0_a0(a,x,y)
                0427           integer(kind=w2f__i4), intent(in) :: a
                0428           type(active), intent(in) :: x
                0429           type(active), intent(inout) :: y
                0430           VECTOR_LOOP_VAR
                0431           VECTOR_LOOP_BEGIN
                0432             y%DELEM=x%DELEM*a
                0433           VECTOR_LOOP_END
                0434         end subroutine
                0435         subroutine sax_d0_a0_a1(a,x,y)
                0436           real(w2f__8), intent(in) :: a
                0437           type(active), intent(in) :: x
                0438           type(active), dimension(:), intent(inout) :: y
                0439           VECTOR_LOOP_VAR
                0440           VECTOR_LOOP_BEGIN
                0441             y%DELEM=x%DELEM*a
                0442           VECTOR_LOOP_END
                0443         end subroutine
                0444         subroutine sax_d0_a0_a2(a,x,y)
                0445           real(w2f__8), intent(in) :: a
                0446           type(active), intent(in) :: x
                0447           type(active), dimension(:,:), intent(inout) :: y
                0448           VECTOR_LOOP_VAR
                0449           VECTOR_LOOP_BEGIN
                0450             y%DELEM = x%DELEM*a
                0451           VECTOR_LOOP_END
                0452         end subroutine
                0453         subroutine sax_d0_a1_a1(a,x,y)
                0454           real(w2f__8), intent(in) :: a
                0455           type(active), dimension(:), intent(in) :: x
                0456           type(active), dimension(:), intent(inout) :: y
                0457           VECTOR_LOOP_VAR
                0458           VECTOR_LOOP_BEGIN
                0459             y%DELEM=x%DELEM*a
                0460           VECTOR_LOOP_END
                0461         end subroutine
                0462         subroutine sax_d0_a2_a2(a,x,y)
                0463           real(w2f__8), intent(in) :: a
                0464           type(active), dimension(:,:), intent(in) :: x
                0465           type(active), dimension(:,:), intent(inout) :: y
                0466           VECTOR_LOOP_VAR
                0467           VECTOR_LOOP_BEGIN
                0468             y%DELEM=x%DELEM*a
                0469           VECTOR_LOOP_END
                0470         end subroutine
                0471         subroutine sax_d0_a3_a3(a,x,y)
                0472           real(w2f__8), intent(in) :: a
                0473           type(active), dimension(:,:,:), intent(in) :: x
                0474           type(active), dimension(:,:,:), intent(inout) :: y
                0475           VECTOR_LOOP_VAR
                0476           VECTOR_LOOP_BEGIN
                0477             y%DELEM=x%DELEM*a
                0478           VECTOR_LOOP_END
                0479         end subroutine
                0480         subroutine sax_d1_a0_a1(a,x,y)
                0481           real(w2f__8), dimension(:), intent(in) :: a
                0482           type(active), intent(in) :: x
                0483           type(active), dimension(:), intent(inout) :: y
                0484           VECTOR_LOOP_VAR
                0485           VECTOR_LOOP_BEGIN
                0486             y%DELEM=x%DELEM*a
                0487           VECTOR_LOOP_END
                0488         end subroutine
                0489         subroutine sax_d1_a1_a1(a,x,y)
                0490           real(w2f__8), dimension(:), intent(in) :: a
                0491           type(active), dimension(:), intent(in) :: x
                0492           type(active), dimension(:), intent(inout) :: y
                0493           VECTOR_LOOP_VAR
                0494           VECTOR_LOOP_BEGIN
                0495             y%DELEM=x%DELEM*a
                0496           VECTOR_LOOP_END
                0497         end subroutine
                0498         subroutine sax_l1_a1_a1(a,x,y)
                0499           integer(kind=w2f__i8), dimension(:), intent(in) :: a
                0500           type(active), dimension(:), intent(in) :: x
                0501           type(active), dimension(:), intent(inout) :: y
                0502           VECTOR_LOOP_VAR
                0503           VECTOR_LOOP_BEGIN
                0504             y%DELEM=x%DELEM*a
                0505           VECTOR_LOOP_END
                0506         end subroutine
                0507         subroutine sax_i1_a1_a1(a,x,y)
                0508           integer(kind=w2f__i4), dimension(:), intent(in) :: a
                0509           type(active), dimension(:), intent(in) :: x
                0510           type(active), dimension(:), intent(inout) :: y
                0511           VECTOR_LOOP_VAR
                0512           VECTOR_LOOP_BEGIN
                0513             y%DELEM=x%DELEM*a
                0514           VECTOR_LOOP_END
                0515         end subroutine
                0516         subroutine sax_d2_a0_a2(a,x,y)
                0517           real(w2f__8), dimension(:,:), intent(in) :: a
                0518           type(active), intent(in) :: x
                0519           type(active), dimension(:,:), intent(inout) :: y
                0520           VECTOR_LOOP_VAR
                0521           VECTOR_LOOP_BEGIN
                0522             y%DELEM=x%DELEM*a
                0523           VECTOR_LOOP_END
                0524         end subroutine
                0525         subroutine sax_d2_a2_a2(a,x,y)
                0526           real(w2f__8), dimension(:,:), intent(in) :: a
                0527           type(active), dimension(:,:), intent(in) :: x
                0528           type(active), dimension(:,:), intent(inout) :: y
                0529           VECTOR_LOOP_VAR
                0530           VECTOR_LOOP_BEGIN
                0531             y%DELEM=x%DELEM*a
                0532           VECTOR_LOOP_END
                0533         end subroutine
                0534 # ifndef DEFAULT_R8
                0535         subroutine sax_r0_a0_a0(a,x,y)
                0536           real(w2f__4), intent(in) :: a
                0537           type(active), intent(in) :: x
                0538           type(active), intent(inout) :: y
                0539           VECTOR_LOOP_VAR
                0540           VECTOR_LOOP_BEGIN
                0541             y%DELEM=x%DELEM*a
                0542           VECTOR_LOOP_END
                0543         end subroutine
                0544 # endif
                0545         !
                0546         ! set derivative of y to be equal to derivative of x
                0547         ! note: making y inout allows for already existing active
                0548         ! variables to become the target of a derivative assignment
                0549         !
                0550         subroutine setderiv_a0_a0(y,x)
                0551           type(active), intent(inout) :: y
                0552           type(active), intent(in) :: x
                0553           VECTOR_LOOP_VAR
                0554           VECTOR_LOOP_BEGIN
                0555             y%DELEM = x%DELEM
                0556           VECTOR_LOOP_END
                0557         end subroutine
                0558         subroutine setderiv_a1_a0(y,x)
                0559           type(active), intent(inout), dimension(:) :: y
                0560           type(active), intent(in) :: x
                0561           VECTOR_LOOP_VAR
                0562           VECTOR_LOOP_BEGIN
                0563             y%DELEM = x%DELEM
                0564           VECTOR_LOOP_END
                0565         end subroutine
                0566         subroutine setderiv_a1_a1(y,x)
                0567           type(active), intent(inout), dimension(:) :: y
                0568           type(active), intent(in), dimension(:) :: x
                0569           VECTOR_LOOP_VAR
                0570           VECTOR_LOOP_BEGIN
                0571             y%DELEM = x%DELEM
                0572           VECTOR_LOOP_END
                0573         end subroutine
                0574         subroutine setderiv_a2_a0(y,x)
                0575           type(active), intent(inout), dimension(:,:) :: y
                0576           type(active), intent(in) :: x
                0577           VECTOR_LOOP_VAR
                0578           VECTOR_LOOP_BEGIN
                0579             y%DELEM = x%DELEM
                0580           VECTOR_LOOP_END
                0581         end subroutine
                0582         subroutine setderiv_a2_a2(y,x)
                0583           type(active), intent(inout), dimension(:,:) :: y
                0584           type(active), intent(in), dimension(:,:) :: x
                0585           VECTOR_LOOP_VAR
                0586           VECTOR_LOOP_BEGIN
                0587             y%DELEM = x%DELEM
                0588           VECTOR_LOOP_END
                0589         end subroutine
                0590         subroutine setderiv_a3_a3(y,x)
                0591           type(active), intent(inout), dimension(:,:,:) :: y
                0592           type(active), intent(in), dimension(:,:,:) :: x
                0593           VECTOR_LOOP_VAR
                0594           VECTOR_LOOP_BEGIN
                0595             y%DELEM = x%DELEM
                0596           VECTOR_LOOP_END
                0597         end subroutine
                0598         !
                0599         ! set the derivative of y to be the negated derivative of x
                0600         ! note: making y inout allows for already existing active
                0601         ! variables to become the target of a derivative assignment
                0602         !
                0603         
                0604         subroutine set_neg_deriv_a0_a0(y,x)
                0605           type(active), intent(inout) :: y
                0606           type(active), intent(in) :: x
                0607           VECTOR_LOOP_VAR
                0608           VECTOR_LOOP_BEGIN
                0609             y%DELEM = -x%DELEM
                0610           VECTOR_LOOP_END
                0611         end subroutine
                0612 
                0613         subroutine set_neg_deriv_a1_a1(y,x)
                0614           type(active), intent(inout), dimension(:) :: y
                0615           type(active), intent(in), dimension(:) :: x
                0616           VECTOR_LOOP_VAR
                0617           VECTOR_LOOP_BEGIN
                0618             y%DELEM = -x%DELEM
                0619           VECTOR_LOOP_END
                0620         end subroutine
                0621 
                0622         !
                0623         ! increment the derivative of y by the derivative of x
                0624         ! note: making y inout allows for already existing active
                0625         ! variables to become the target of a derivative assignment
                0626         !
                0627         
                0628         subroutine inc_deriv_a0_a0(y,x)
                0629           type(active), intent(inout) :: y
                0630           type(active), intent(in) :: x
                0631           VECTOR_LOOP_VAR
                0632           VECTOR_LOOP_BEGIN
                0633             y%DELEM=y%DELEM + x%DELEM
                0634           VECTOR_LOOP_END
                0635         end subroutine
                0636 
                0637         subroutine inc_deriv_a1_a1(y,x)
                0638           type(active), intent(inout), dimension(:) :: y
                0639           type(active), intent(in), dimension(:) :: x
                0640           VECTOR_LOOP_VAR
                0641           VECTOR_LOOP_BEGIN
                0642             y%DELEM=y%DELEM + x%DELEM
                0643           VECTOR_LOOP_END
                0644         end subroutine
                0645 
                0646         subroutine inc_deriv_a2_a2(y,x)
                0647           type(active), intent(inout), dimension(:,:) :: y
                0648           type(active), intent(in), dimension(:,:) :: x
                0649           VECTOR_LOOP_VAR
                0650           VECTOR_LOOP_BEGIN
                0651             y%DELEM=y%DELEM + x%DELEM
                0652           VECTOR_LOOP_END
                0653         end subroutine
                0654 
                0655         !
                0656         ! decrement the derivative of y by the derivative of x
                0657         ! note: making y inout allows for already existing active
                0658         ! variables to become the target of a derivative assignment
                0659         !
                0660         
                0661         subroutine dec_deriv_a0_a0(y,x)
                0662           type(active), intent(inout) :: y
                0663           type(active), intent(in) :: x
                0664           VECTOR_LOOP_VAR
                0665           VECTOR_LOOP_BEGIN
                0666             y%DELEM=y%DELEM - x%DELEM
                0667           VECTOR_LOOP_END
                0668         end subroutine
                0669 
                0670         subroutine dec_deriv_a1_a1(y,x)
                0671           type(active), intent(inout), dimension(:) :: y
                0672           type(active), intent(in), dimension(:) :: x
                0673           VECTOR_LOOP_VAR
                0674           VECTOR_LOOP_BEGIN
                0675             y%DELEM=y%DELEM - x%DELEM
                0676           VECTOR_LOOP_END
                0677         end subroutine
                0678         
                0679         subroutine dec_deriv_a2_a2(y,x)
                0680           type(active), intent(inout), dimension(:,:) :: y
                0681           type(active), intent(in), dimension(:,:) :: x
                0682           VECTOR_LOOP_VAR
                0683           VECTOR_LOOP_BEGIN
                0684             y%DELEM=y%DELEM - x%DELEM
                0685           VECTOR_LOOP_END
                0686         end subroutine
                0687         
                0688         !
                0689         ! set derivative components to 0.0
                0690         !
                0691         subroutine zero_deriv_a0(x)
                0692           type(active), intent(inout) :: x
                0693           VECTOR_LOOP_VAR
                0694           VECTOR_LOOP_BEGIN
                0695              x%DELEM=0.0d0
                0696           VECTOR_LOOP_END 
                0697         end subroutine
                0698 
                0699         subroutine zero_deriv_a1(x)
                0700           type(active), dimension(:), intent(inout) :: x
                0701           VECTOR_LOOP_VAR
                0702           VECTOR_LOOP_BEGIN
                0703              x%DELEM=0.0d0
                0704           VECTOR_LOOP_END 
                0705         end subroutine
                0706 
                0707         subroutine zero_deriv_a2(x)
                0708           type(active), dimension(:,:), intent(inout) :: x
                0709           VECTOR_LOOP_VAR
                0710           VECTOR_LOOP_BEGIN
                0711              x%DELEM=0.0d0
                0712           VECTOR_LOOP_END 
                0713         end subroutine
                0714 
                0715         subroutine zero_deriv_a3(x)
                0716           type(active), dimension(:,:,:), intent(inout) :: x
                0717           VECTOR_LOOP_VAR
                0718           VECTOR_LOOP_BEGIN
                0719              x%DELEM=0.0d0
                0720           VECTOR_LOOP_END 
                0721         end subroutine
                0722 
                0723         subroutine zero_deriv_a4(x)
                0724           type(active), dimension(:,:,:,:), intent(inout) :: x
                0725           VECTOR_LOOP_VAR
                0726           VECTOR_LOOP_BEGIN
                0727              x%DELEM=0.0d0
                0728           VECTOR_LOOP_END 
                0729         end subroutine
                0730 
                0731 #endif
                0732         !
                0733         ! conversions
                0734         !
                0735         subroutine convert_d0_a0(convertTo, convertFrom)
                0736           real(w2f__8), intent(out) :: convertTo
                0737           type(active), intent(in) :: convertFrom
                0738           convertTo=convertFrom%v
                0739         end subroutine
                0740         subroutine convert_d1_a1(convertTo, convertFrom)
                0741           real(w2f__8), dimension(:), intent(out) :: convertTo
                0742           type(active), dimension(:), intent(in) :: convertFrom
                0743           convertTo=convertFrom%v
                0744         end subroutine
                0745         subroutine convert_d2_a2(convertTo, convertFrom)
                0746           real(w2f__8), dimension(:,:), intent(out) :: convertTo
                0747           type(active), dimension(:,:), intent(in) :: convertFrom
                0748           convertTo=convertFrom%v
                0749         end subroutine
                0750         subroutine convert_d3_a3(convertTo, convertFrom)
                0751           real(w2f__8), dimension(:,:,:), intent(out) :: convertTo
                0752           type(active), dimension(:,:,:), intent(in) :: convertFrom
                0753           convertTo=convertFrom%v
                0754         end subroutine
                0755         subroutine convert_d4_a4(convertTo, convertFrom)
                0756           real(w2f__8), dimension(:,:,:,:), intent(out) :: convertTo
                0757           type(active), dimension(:,:,:,:), intent(in) :: convertFrom
                0758           convertTo=convertFrom%v
                0759         end subroutine
                0760         subroutine convert_d5_a5(convertTo, convertFrom)
                0761           real(w2f__8), dimension(:,:,:,:,:), intent(out) :: convertTo
                0762           type(active), dimension(:,:,:,:,:), intent(in) :: convertFrom
                0763           convertTo=convertFrom%v
                0764         end subroutine
                0765         subroutine convert_d6_a6(convertTo, convertFrom)
                0766           real(w2f__8), dimension(:,:,:,:,:,:), intent(out) :: convertTo
                0767           type(active), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
                0768           convertTo=convertFrom%v
                0769         end subroutine
                0770         subroutine convert_d7_a7(convertTo, convertFrom)
                0771           real(w2f__8), dimension(:,:,:,:,:,:,:), intent(out) :: convertTo
                0772           type(active), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
                0773           convertTo=convertFrom%v
                0774         end subroutine
                0775 
                0776         subroutine convert_a0_d0(convertTo, convertFrom)
                0777           type(active), intent(inout) :: convertTo
                0778           real(w2f__8), intent(in) :: convertFrom
                0779           convertTo%v=convertFrom
                0780         end subroutine
                0781         subroutine convert_a1_d1(convertTo, convertFrom)
                0782           type(active), dimension(:), intent(inout) :: convertTo
                0783           real(w2f__8), dimension(:), intent(in) :: convertFrom
                0784           convertTo%v=convertFrom
                0785         end subroutine
                0786         subroutine convert_a2_d2(convertTo, convertFrom)
                0787           type(active), dimension(:,:), intent(inout) :: convertTo
                0788           real(w2f__8), dimension(:,:), intent(in) :: convertFrom
                0789           convertTo%v=convertFrom
                0790         end subroutine
                0791         subroutine convert_a3_d3(convertTo, convertFrom)
                0792           type(active), dimension(:,:,:), intent(inout) :: convertTo
                0793           real(w2f__8), dimension(:,:,:), intent(in) :: convertFrom
                0794           convertTo%v=convertFrom
                0795         end subroutine
                0796         subroutine convert_a4_d4(convertTo, convertFrom)
                0797           type(active), dimension(:,:,:,:), intent(inout) :: convertTo
                0798           real(w2f__8), dimension(:,:,:,:), intent(in) :: convertFrom
                0799           convertTo%v=convertFrom
                0800         end subroutine
                0801         subroutine convert_a5_d5(convertTo, convertFrom)
                0802           type(active), dimension(:,:,:,:,:), intent(inout) :: convertTo
                0803           real(w2f__8), dimension(:,:,:,:,:), intent(in) :: convertFrom
                0804           convertTo%v=convertFrom
                0805         end subroutine
                0806         subroutine convert_a6_d6(convertTo, convertFrom)
                0807           type(active), dimension(:,:,:,:,:,:), intent(inout) :: convertTo
                0808           real(w2f__8), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
                0809           convertTo%v=convertFrom
                0810         end subroutine
                0811         subroutine convert_a7_d7(convertTo, convertFrom)
                0812           type(active), dimension(:,:,:,:,:,:,:), intent(inout) :: convertTo
                0813           real(w2f__8), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
                0814           convertTo%v=convertFrom
                0815         end subroutine
                0816 #ifndef DEFAULT_R8
                0817         subroutine convert_r0_a0(convertTo, convertFrom)
                0818           real(w2f__4), intent(out) :: convertTo
                0819           type(active), intent(in) :: convertFrom
                0820           convertTo=convertFrom%v
                0821         end subroutine
                0822         subroutine convert_r1_a1(convertTo, convertFrom)
                0823           real(w2f__4), dimension(:), intent(out) :: convertTo
                0824           type(active), dimension(:), intent(in) :: convertFrom
                0825           convertTo=convertFrom%v
                0826         end subroutine
                0827         subroutine convert_r2_a2(convertTo, convertFrom)
                0828           real(w2f__4), dimension(:,:), intent(out) :: convertTo
                0829           type(active), dimension(:,:), intent(in) :: convertFrom
                0830           convertTo=convertFrom%v
                0831         end subroutine
                0832         subroutine convert_r3_a3(convertTo, convertFrom)
                0833           real(w2f__4), dimension(:,:,:), intent(out) :: convertTo
                0834           type(active), dimension(:,:,:), intent(in) :: convertFrom
                0835           convertTo=convertFrom%v
                0836         end subroutine
                0837         subroutine convert_r4_a4(convertTo, convertFrom)
                0838           real(w2f__4), dimension(:,:,:,:), intent(out) :: convertTo
                0839           type(active), dimension(:,:,:,:), intent(in) :: convertFrom
                0840           convertTo=convertFrom%v
                0841         end subroutine
                0842         subroutine convert_r5_a5(convertTo, convertFrom)
                0843           real(w2f__4), dimension(:,:,:,:,:), intent(out) :: convertTo
                0844           type(active), dimension(:,:,:,:,:), intent(in) :: convertFrom
                0845           convertTo=convertFrom%v
                0846         end subroutine
                0847         subroutine convert_r6_a6(convertTo, convertFrom)
                0848           real(w2f__4), dimension(:,:,:,:,:,:), intent(out) :: convertTo
                0849           type(active), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
                0850           convertTo=convertFrom%v
                0851         end subroutine
                0852         subroutine convert_r7_a7(convertTo, convertFrom)
                0853           real(w2f__4), dimension(:,:,:,:,:,:,:), intent(out) :: convertTo
                0854           type(active), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
                0855           convertTo=convertFrom%v
                0856         end subroutine
                0857 
                0858         subroutine convert_a0_r0(convertTo, convertFrom)
                0859           type(active), intent(inout) :: convertTo
                0860           real(w2f__4), intent(in) :: convertFrom
                0861           convertTo%v=convertFrom
                0862         end subroutine
                0863         subroutine convert_a1_r1(convertTo, convertFrom)
                0864           type(active), dimension(:), intent(inout) :: convertTo
                0865           real(w2f__4), dimension(:), intent(in) :: convertFrom
                0866           convertTo%v=convertFrom
                0867         end subroutine
                0868         subroutine convert_a2_r2(convertTo, convertFrom)
                0869           type(active), dimension(:,:), intent(inout) :: convertTo
                0870           real(w2f__4), dimension(:,:), intent(in) :: convertFrom
                0871           convertTo%v=convertFrom
                0872         end subroutine
                0873         subroutine convert_a3_r3(convertTo, convertFrom)
                0874           type(active), dimension(:,:,:), intent(inout) :: convertTo
                0875           real(w2f__4), dimension(:,:,:), intent(in) :: convertFrom
                0876           convertTo%v=convertFrom
                0877         end subroutine
                0878         subroutine convert_a4_r4(convertTo, convertFrom)
                0879           type(active), dimension(:,:,:,:), intent(inout) :: convertTo
                0880           real(w2f__4), dimension(:,:,:,:), intent(in) :: convertFrom
                0881           convertTo%v=convertFrom
                0882         end subroutine
                0883         subroutine convert_a5_r5(convertTo, convertFrom)
                0884           type(active), dimension(:,:,:,:,:), intent(inout) :: convertTo
                0885           real(w2f__4), dimension(:,:,:,:,:), intent(in) :: convertFrom
                0886           convertTo%v=convertFrom
                0887         end subroutine
                0888         subroutine convert_a6_r6(convertTo, convertFrom)
                0889           type(active), dimension(:,:,:,:,:,:), intent(inout) :: convertTo
                0890           real(w2f__4), dimension(:,:,:,:,:,:), intent(in) :: convertFrom
                0891           convertTo%v=convertFrom
                0892         end subroutine
                0893         subroutine convert_a7_r7(convertTo, convertFrom)
                0894           type(active), dimension(:,:,:,:,:,:,:), intent(inout) :: convertTo
                0895           real(w2f__4), dimension(:,:,:,:,:,:,:), intent(in) :: convertFrom
                0896           convertTo%v=convertFrom
                0897         end subroutine
                0898 #endif
                0899         !
                0900         ! allocations
                0901         !
                0902         subroutine allocateMatching_d1_d1(toBeAllocated,allocateMatching)
                0903           implicit none
                0904           real(w2f__8), dimension(:), allocatable :: toBeAllocated
                0905           real(w2f__8), dimension(:) :: allocateMatching
                0906           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0907           allocate(toBeAllocated(size(allocateMatching)))
                0908         end subroutine
                0909         subroutine allocateMatching_a1_d1(toBeAllocated,allocateMatching)
                0910           implicit none
                0911           type(active), dimension(:), allocatable :: toBeAllocated
                0912           real(w2f__8), dimension(:) :: allocateMatching
                0913           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0914           allocate(toBeAllocated(size(allocateMatching)))
                0915         end subroutine
                0916         subroutine allocateMatching_d1_a1(toBeAllocated,allocateMatching)
                0917           implicit none
                0918           real(w2f__8), dimension(:), allocatable :: toBeAllocated
                0919           type(active), dimension(:) :: allocateMatching
                0920           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0921           allocate(toBeAllocated(size(allocateMatching)))
                0922         end subroutine
                0923         subroutine allocateMatching_a1_a1(toBeAllocated,allocateMatching)
                0924           implicit none
                0925           type(active), dimension(:), allocatable :: toBeAllocated
                0926           type(active), dimension(:) :: allocateMatching
                0927           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0928           allocate(toBeAllocated(size(allocateMatching)))
                0929         end subroutine
                0930         subroutine allocateMatching_d2_d2(toBeAllocated,allocateMatching)
                0931           implicit none
                0932           real(w2f__8), dimension(:,:), allocatable :: toBeAllocated
                0933           real(w2f__8), dimension(:,:) :: allocateMatching
                0934           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0935           allocate(toBeAllocated(size(allocateMatching,1), &
                0936                size(allocateMatching,2)))
                0937         end subroutine
                0938         subroutine allocateMatching_a2_d2(toBeAllocated,allocateMatching)
                0939           implicit none
                0940           type(active), dimension(:,:), allocatable :: toBeAllocated
                0941           real(w2f__8), dimension(:,:) :: allocateMatching
                0942           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0943           allocate(toBeAllocated(size(allocateMatching,1), &
                0944                size(allocateMatching,2)))
                0945         end subroutine
                0946         subroutine allocateMatching_d2_a2(toBeAllocated,allocateMatching)
                0947           implicit none
                0948           real(w2f__8), dimension(:,:), allocatable :: toBeAllocated
                0949           type(active), dimension(:,:) :: allocateMatching
                0950           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0951           allocate(toBeAllocated(size(allocateMatching,1), &
                0952                size(allocateMatching,2)))
                0953         end subroutine
                0954         subroutine allocateMatching_a2_a2(toBeAllocated,allocateMatching)
                0955           implicit none
                0956           type(active), dimension(:,:), allocatable :: toBeAllocated
                0957           type(active), dimension(:,:) :: allocateMatching
                0958           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0959           allocate(toBeAllocated(size(allocateMatching,1), &
                0960                size(allocateMatching,2)))
                0961         end subroutine
                0962         subroutine allocateMatching_d3_d3(toBeAllocated,allocateMatching)
                0963           implicit none
                0964           real(w2f__8), dimension(:,:,:), allocatable :: toBeAllocated
                0965           real(w2f__8), dimension(:,:,:) :: allocateMatching
                0966           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0967           allocate(toBeAllocated(size(allocateMatching,1), &
                0968                size(allocateMatching,2),&
                0969                size(allocateMatching,3)))
                0970         end subroutine
                0971         subroutine allocateMatching_a3_d3(toBeAllocated,allocateMatching)
                0972           implicit none
                0973           type(active), dimension(:,:,:), allocatable :: toBeAllocated
                0974           real(w2f__8), dimension(:,:,:) :: allocateMatching
                0975           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0976           allocate(toBeAllocated(size(allocateMatching,1), &
                0977                size(allocateMatching,2),&
                0978                size(allocateMatching,3)))
                0979         end subroutine
                0980         subroutine allocateMatching_d3_a3(toBeAllocated,allocateMatching)
                0981           implicit none
                0982           real(w2f__8), dimension(:,:,:), allocatable :: toBeAllocated
                0983           type(active), dimension(:,:,:) :: allocateMatching
                0984           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0985           allocate(toBeAllocated(size(allocateMatching,1), &
                0986                size(allocateMatching,2),&
                0987                size(allocateMatching,3)))
                0988         end subroutine
                0989         subroutine allocateMatching_a3_a3(toBeAllocated,allocateMatching)
                0990           implicit none
                0991           type(active), dimension(:,:,:), allocatable :: toBeAllocated
                0992           type(active), dimension(:,:,:) :: allocateMatching
                0993           if (allocated(toBeAllocated)) deallocate(toBeAllocated)
                0994           allocate(toBeAllocated(size(allocateMatching,1), &
                0995                size(allocateMatching,2),&
                0996                size(allocateMatching,3)))
                0997         end subroutine
                0998         subroutine allocateMatching_a4_a4(toBeAllocated,allocateMatching)
                0999           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