Back to home page

MITgcm

 
 

    


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

view on githubraw file Latest commit 279dc77b on 2015-07-03 21:33:55 UTC
8702af1f36 Patr*0001 C taping --------------------------------------------
                0002 
                0003 
                0004       subroutine push_s0(x)
                0005 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0006       use OAD_tape
8702af1f36 Patr*0007       implicit none
                0008       double precision :: x
                0009 C $OpenAD$ END DECLS
                0010       if(oad_dt_sz .lt. oad_dt_ptr) call oad_dt_grow()
                0011       oad_dt(oad_dt_ptr)=x; oad_dt_ptr=oad_dt_ptr+1
                0012       end subroutine 
                0013 
                0014       subroutine pop_s0(x)
                0015 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0016       use OAD_tape
8702af1f36 Patr*0017       implicit none
                0018       double precision :: x
                0019 C $OpenAD$ END DECLS
                0020       oad_dt_ptr=oad_dt_ptr-1
                0021       x=oad_dt(oad_dt_ptr)
                0022       end subroutine
                0023 
9c5e597b65 Jean*0024       subroutine push_s1(x)
                0025 C $OpenAD$ INLINE DECLS
                0026       use OAD_tape
                0027       implicit none
                0028       double precision :: x(:)
                0029 C $OpenAD$ END DECLS
                0030       oad_chunk_size=size(x,1)
                0031       if(oad_dt_sz .lt. oad_dt_ptr+oad_chunk_size)
                0032      + call oad_dt_grow()
                0033       oad_dt(oad_dt_ptr:oad_dt_ptr+oad_chunk_size-1)=
                0034      +x
                0035       oad_dt_ptr=oad_dt_ptr+oad_chunk_size
                0036       end subroutine 
                0037 
                0038       subroutine pop_s1(x)
                0039 C $OpenAD$ INLINE DECLS
                0040       use OAD_tape
                0041       implicit none
                0042       double precision :: x(:)
                0043 C $OpenAD$ END DECLS
                0044       oad_chunk_size=size(x,1)
                0045       oad_dt_ptr=oad_dt_ptr-oad_chunk_size
                0046       x=oad_dt(oad_dt_ptr:oad_dt_ptr+oad_chunk_size-1)
                0047       end subroutine
                0048 
                0049       subroutine push_s2(x)
                0050 C $OpenAD$ INLINE DECLS
                0051       use OAD_tape
                0052       implicit none
                0053       double precision :: x(:,:)
                0054 C $OpenAD$ END DECLS
                0055       oad_chunk_size=size(x,1)*size(x,2)
                0056       if(oad_dt_sz .lt. oad_dt_ptr+oad_chunk_size) 
                0057      + call oad_dt_grow()
                0058       oad_dt(oad_dt_ptr:oad_dt_ptr+oad_chunk_size-1)=
                0059      +reshape(x,(/oad_chunk_size/))
                0060       oad_dt_ptr=oad_dt_ptr+oad_chunk_size
                0061       end subroutine 
                0062 
                0063       subroutine pop_s2(x)
                0064 C $OpenAD$ INLINE DECLS
                0065       use OAD_tape
                0066       implicit none
                0067       double precision :: x(:,:)
                0068 C $OpenAD$ END DECLS
                0069       oad_chunk_size=size(x,1)*size(x,2)
                0070       oad_dt_ptr=oad_dt_ptr-oad_chunk_size
                0071         x=reshape(oad_dt(oad_dt_ptr:oad_dt_ptr+oad_chunk_size-1),
                0072      +shape(x))
                0073       end subroutine
                0074 
8702af1f36 Patr*0075       subroutine apush(x)
                0076 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0077       use OAD_tape
                0078       use OAD_active
8702af1f36 Patr*0079       implicit none
                0080       type(active) :: x
                0081 C $OpenAD$ END DECLS
                0082       if(oad_dt_sz .lt. oad_dt_ptr) call oad_dt_grow()
                0083       oad_dt(oad_dt_ptr)=x%v; oad_dt_ptr=oad_dt_ptr+1
                0084       end subroutine 
                0085 
                0086       subroutine apop(x)
                0087 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0088       use OAD_tape
                0089       use OAD_active
8702af1f36 Patr*0090       implicit none
                0091       type(active) :: x
                0092 C $OpenAD$ END DECLS
                0093       oad_dt_ptr=oad_dt_ptr-1
                0094       x%v=oad_dt(oad_dt_ptr)
                0095       end subroutine
                0096 
                0097       subroutine push_i_s0(x)
                0098 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0099       use OAD_tape
8702af1f36 Patr*0100       implicit none
                0101       integer :: x
                0102 C $OpenAD$ END DECLS
                0103       if(oad_it_sz .lt. oad_it_ptr) call oad_it_grow()
                0104       oad_it(oad_it_ptr)=x; oad_it_ptr=oad_it_ptr+1
                0105       end subroutine 
                0106 
                0107       subroutine pop_i_s0(x)
                0108 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0109       use OAD_tape
8702af1f36 Patr*0110       implicit none
                0111       integer :: x
                0112 C $OpenAD$ END DECLS
                0113       oad_it_ptr=oad_it_ptr-1
                0114       x=oad_it(oad_it_ptr)
                0115       end subroutine
                0116 
9c5e597b65 Jean*0117       subroutine push_i_s1(x)
                0118 C $OpenAD$ INLINE DECLS
                0119       use OAD_tape
                0120       implicit none
                0121       integer :: x(:)
                0122 C $OpenAD$ END DECLS
                0123       oad_chunk_size=size(x,1)
                0124       if(oad_it_sz .lt. oad_it_ptr+oad_chunk_size) 
                0125      +call oad_it_grow()
                0126       oad_it(oad_it_ptr:oad_it_ptr+oad_chunk_size-1)=
                0127      +x 
                0128       oad_it_ptr=oad_it_ptr+oad_chunk_size
                0129       end subroutine 
                0130 
                0131       subroutine pop_i_s1(x)
                0132 C $OpenAD$ INLINE DECLS
                0133       use OAD_tape
                0134       implicit none
                0135       integer :: x(:)
                0136 C $OpenAD$ END DECLS
                0137       oad_chunk_size=size(x,1)
                0138       oad_it_ptr=oad_it_ptr-oad_chunk_size
                0139       x=oad_it(oad_it_ptr:oad_it_ptr+oad_chunk_size-1)
                0140       end subroutine
                0141 
                0142       subroutine push_i_s2(x)
                0143 C $OpenAD$ INLINE DECLS
                0144       use OAD_tape
                0145       implicit none
                0146       integer :: x(:,:)
                0147 C $OpenAD$ END DECLS
                0148       oad_chunk_size=size(x,1)*size(x,2)
                0149       if(oad_it_sz .lt. oad_it_ptr+oad_chunk_size) 
                0150      + call oad_it_grow()
                0151       oad_it(oad_it_ptr:oad_it_ptr+oad_chunk_size-1)=
                0152      +reshape(x,(/oad_chunk_size/))
                0153       oad_it_ptr=oad_it_ptr+oad_chunk_size
                0154       end subroutine 
                0155 
                0156       subroutine pop_i_s2(x)
                0157 C $OpenAD$ INLINE DECLS
                0158       use OAD_tape
                0159       implicit none
                0160       integer :: x(:,:)
                0161 C $OpenAD$ END DECLS
                0162       oad_chunk_size=size(x,1)*size(x,2)
                0163       oad_it_ptr=oad_it_ptr-oad_chunk_size
                0164         x=reshape(oad_it(oad_it_ptr:oad_it_ptr+oad_chunk_size-1),
                0165      +shape(x))
                0166       end subroutine
                0167 
8702af1f36 Patr*0168       subroutine push_b(x)
                0169 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0170       use OAD_tape
8702af1f36 Patr*0171       implicit none
                0172       logical :: x
                0173 C $OpenAD$ END DECLS
                0174       if(oad_lt_sz .lt. oad_lt_ptr) call oad_lt_grow()
                0175       oad_lt(oad_lt_ptr)=x; oad_lt_ptr=oad_lt_ptr+1
                0176       end subroutine 
                0177 
                0178       subroutine pop_b(x)
                0179 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0180       use OAD_tape
8702af1f36 Patr*0181       implicit none
                0182       logical :: x
                0183 C $OpenAD$ END DECLS
                0184       oad_lt_ptr=oad_lt_ptr-1
                0185       x=oad_lt(oad_lt_ptr)
                0186       end subroutine
                0187 
7d8fa04723 Jean*0188       subroutine push_s(s)
8702af1f36 Patr*0189 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0190       use OAD_tape
8702af1f36 Patr*0191       implicit none
7d8fa04723 Jean*0192       character*(80) :: s
8702af1f36 Patr*0193 C $OpenAD$ END DECLS
                0194       if(oad_st_sz .lt. oad_st_ptr) call oad_st_grow()
7d8fa04723 Jean*0195       oad_st(oad_st_ptr)=s; oad_st_ptr=oad_st_ptr+1
8702af1f36 Patr*0196       end subroutine 
                0197 
7d8fa04723 Jean*0198       subroutine pop_s(s)
8702af1f36 Patr*0199 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0200       use OAD_tape
8702af1f36 Patr*0201       implicit none
7d8fa04723 Jean*0202       character*(80) :: s
8702af1f36 Patr*0203 C $OpenAD$ END DECLS
                0204       oad_st_ptr=oad_st_ptr-1
7d8fa04723 Jean*0205       s=oad_st(oad_st_ptr)
8702af1f36 Patr*0206       end subroutine
9c5e597b65 Jean*0207 
                0208 C ----------------------- Propagation -----------------------
                0209 
8702af1f36 Patr*0210       subroutine saxpy(a,x,y)
                0211 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0212       use OAD_active
7d8fa04723 Jean*0213       implicit none
8702af1f36 Patr*0214       double precision, intent(in) :: a
                0215       type(active), intent(in) :: x
                0216       type(active), intent(inout) :: y
                0217 C $OpenAD$ END DECLS
                0218       y%d=y%d+x%d*(a)
                0219       end subroutine
                0220 
                0221       subroutine zeroderiv(x)
                0222 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0223       use OAD_active
7d8fa04723 Jean*0224       implicit none
8702af1f36 Patr*0225       type(active), intent(out) :: x
                0226 C $OpenAD$ END DECLS
                0227       x%d=0.0d0
                0228       end subroutine
                0229 
                0230       subroutine setderiv(y,x)
                0231 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0232       use OAD_active
7d8fa04723 Jean*0233       implicit none
8702af1f36 Patr*0234       type(active), intent(out) :: x
                0235       type(active), intent(in) :: y
                0236 C $OpenAD$ END DECLS
                0237       x%d=y%d
                0238       end subroutine
                0239 
                0240       subroutine incderiv(y,x)
                0241 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0242       use OAD_active
7d8fa04723 Jean*0243       implicit none
8702af1f36 Patr*0244       type(active), intent(out) :: x
                0245       type(active), intent(in) :: y
                0246 C $OpenAD$ END DECLS
                0247       x%d=x%d+y%d
                0248       end subroutine
                0249 
                0250       subroutine decderiv(y,x)
                0251 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0252       use OAD_active
7d8fa04723 Jean*0253       implicit none
8702af1f36 Patr*0254       type(active), intent(out) :: x
                0255       type(active), intent(in) :: y
                0256 C $OpenAD$ END DECLS
                0257       x%d = x%d - y%d
                0258       end subroutine decderiv
                0259 
                0260 C Checkpointing stuff ---------------------------------------
                0261 
                0262 C reals -----------------------------------------------------
                0263       subroutine cp_arg_store_real_scalar(x)
                0264 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0265         use OAD_cp
                0266         implicit none
                0267         double precision :: x
8702af1f36 Patr*0268 C $OpenAD$ END DECLS
                0269 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0270         !write(standardmessageunit,*)'OAD: cp write x ', x
8702af1f36 Patr*0271 #endif
9c5e597b65 Jean*0272         write(unit=cp_io_unit) x
8702af1f36 Patr*0273       end subroutine 
                0274 
                0275       subroutine cp_arg_restore_real_scalar(x)
                0276 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0277         use OAD_cp
                0278         implicit none
                0279         double precision :: x
8702af1f36 Patr*0280 C $OpenAD$ END DECLS
9c5e597b65 Jean*0281         read(unit=cp_io_unit) x
8702af1f36 Patr*0282 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0283         !write(standardmessageunit,*)'OAD: cp read x ', x
8702af1f36 Patr*0284 #endif
                0285       end subroutine 
                0286 
                0287       subroutine cp_arg_store_real_scalar_a(x)
                0288 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0289         use OAD_active
                0290         use OAD_cp
                0291         implicit none
                0292         type(active) :: x
8702af1f36 Patr*0293 C $OpenAD$ END DECLS
                0294 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0295         !write(standardmessageunit,*)'OAD: cp write x ', x%v
8702af1f36 Patr*0296 #endif
9c5e597b65 Jean*0297         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0298       end subroutine 
                0299 
                0300       subroutine cp_arg_restore_real_scalar_a(x)
                0301 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0302         use OAD_active
                0303         use OAD_cp
                0304         implicit none
                0305         type(active) :: x
8702af1f36 Patr*0306 C $OpenAD$ END DECLS
9c5e597b65 Jean*0307         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0308 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0309         !write(standardmessageunit,*)'OAD: cp read x ', x%v
8702af1f36 Patr*0310 #endif
                0311       end subroutine 
7d8fa04723 Jean*0312 
8702af1f36 Patr*0313       subroutine cp_arg_store_real_vector(x)
                0314 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0315         use OAD_cp
                0316         implicit none
                0317         double precision, dimension(:) :: x
8702af1f36 Patr*0318 C $OpenAD$ END DECLS
                0319 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0320         !write(standardmessageunit,*)'OAD: cp write x ', x(1)
8702af1f36 Patr*0321 #endif
9c5e597b65 Jean*0322         write(unit=cp_io_unit) x
8702af1f36 Patr*0323       end subroutine 
                0324 
                0325       subroutine cp_arg_restore_real_vector(x)
                0326 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0327         use OAD_cp
                0328         implicit none
                0329         double precision, dimension(:) :: x
8702af1f36 Patr*0330 C $OpenAD$ END DECLS
9c5e597b65 Jean*0331         read(unit=cp_io_unit) x
8702af1f36 Patr*0332 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0333         !write(standardmessageunit,*)'OAD: cp read x ', x(1)
8702af1f36 Patr*0334 #endif
                0335       end subroutine 
                0336 
                0337       subroutine cp_arg_store_real_vector_a(x)
                0338 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0339         use OAD_active
                0340         use OAD_cp
                0341         implicit none
                0342         type(active), dimension(:) :: x
8702af1f36 Patr*0343 C $OpenAD$ END DECLS
                0344 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0345         !write(standardmessageunit,*)'OAD: cp write x ', x(1)%v
8702af1f36 Patr*0346 #endif
9c5e597b65 Jean*0347         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0348       end subroutine 
                0349 
                0350       subroutine cp_arg_restore_real_vector_a(x)
                0351 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0352         use OAD_active
                0353         use OAD_cp
                0354         implicit none
                0355         type(active), dimension(:) :: x
8702af1f36 Patr*0356 C $OpenAD$ END DECLS
9c5e597b65 Jean*0357         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0358 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0359         !write(standardmessageunit,*)'OAD: cp read x ', x(1)%v
8702af1f36 Patr*0360 #endif
                0361       end subroutine 
                0362 
                0363       subroutine cp_arg_store_real_matrix(x)
                0364 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0365         use OAD_cp
                0366         implicit none
                0367         double precision, dimension(:,:) :: x
8702af1f36 Patr*0368 C $OpenAD$ END DECLS
                0369 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0370         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1)
8702af1f36 Patr*0371 #endif
9c5e597b65 Jean*0372         write(unit=cp_io_unit) x
8702af1f36 Patr*0373       end subroutine 
                0374 
                0375       subroutine cp_arg_restore_real_matrix(x)
                0376 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0377         use OAD_cp
                0378         implicit none
                0379         double precision, dimension(:,:) :: x
8702af1f36 Patr*0380 C $OpenAD$ END DECLS
9c5e597b65 Jean*0381         read(unit=cp_io_unit) x
8702af1f36 Patr*0382 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0383         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1)
8702af1f36 Patr*0384 #endif
                0385       end subroutine 
                0386 
                0387       subroutine cp_arg_store_real_matrix_a(x)
                0388 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0389         use OAD_active
                0390         use OAD_cp
                0391         implicit none
                0392         type(active), dimension(:,:) :: x
8702af1f36 Patr*0393 C $OpenAD$ END DECLS
                0394 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0395         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1)%v
8702af1f36 Patr*0396 #endif
9c5e597b65 Jean*0397         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0398       end subroutine 
                0399 
                0400       subroutine cp_arg_restore_real_matrix_a(x)
                0401 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0402         use OAD_active
                0403         use OAD_cp
                0404         implicit none
                0405         type(active), dimension(:,:) :: x
8702af1f36 Patr*0406 C $OpenAD$ END DECLS
9c5e597b65 Jean*0407         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0408 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0409         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1)%v
8702af1f36 Patr*0410 #endif
                0411       end subroutine 
                0412 
                0413       subroutine cp_arg_store_real_three_tensor(x)
                0414 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0415         use OAD_cp
                0416         implicit none
                0417         double precision, dimension(:,:,:) :: x
8702af1f36 Patr*0418 C $OpenAD$ END DECLS
                0419 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0420         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1)
8702af1f36 Patr*0421 #endif
9c5e597b65 Jean*0422         write(unit=cp_io_unit) x
8702af1f36 Patr*0423       end subroutine 
                0424 
9c5e597b65 Jean*0425       subroutine cp_arg_restore_real_three_tensor(x)
8702af1f36 Patr*0426 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0427         use OAD_cp
                0428         implicit none
                0429         double precision, dimension(:,:,:) :: x
8702af1f36 Patr*0430 C $OpenAD$ END DECLS
9c5e597b65 Jean*0431         read(unit=cp_io_unit) x
8702af1f36 Patr*0432 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0433         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1)
8702af1f36 Patr*0434 #endif
                0435       end subroutine 
                0436 
9c5e597b65 Jean*0437       subroutine cp_arg_store_real_three_tensor_a(x)
8702af1f36 Patr*0438 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0439         use OAD_active
                0440         use OAD_cp
                0441         implicit none
                0442         type(active), dimension(:,:,:) :: x
8702af1f36 Patr*0443 C $OpenAD$ END DECLS
                0444 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0445         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1)%v
8702af1f36 Patr*0446 #endif
9c5e597b65 Jean*0447         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0448       end subroutine 
                0449 
                0450       subroutine cp_arg_restore_real_three_tensor_a(x)
                0451 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0452         use OAD_active
                0453         use OAD_cp
                0454         implicit none
                0455         type(active), dimension(:,:,:) :: x
8702af1f36 Patr*0456 C $OpenAD$ END DECLS
7d8fa04723 Jean*0457         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0458 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0459         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1)%v
8702af1f36 Patr*0460 #endif
                0461       end subroutine 
                0462 
                0463       subroutine cp_arg_store_real_four_tensor(x)
                0464 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0465         use OAD_cp
                0466         implicit none
                0467         double precision, dimension(:,:,:,:) :: x
8702af1f36 Patr*0468 C $OpenAD$ END DECLS
                0469 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0470         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1)
8702af1f36 Patr*0471 #endif
9c5e597b65 Jean*0472         write(unit=cp_io_unit) x
8702af1f36 Patr*0473       end subroutine 
                0474 
9c5e597b65 Jean*0475       subroutine cp_arg_restore_real_four_tensor(x)
8702af1f36 Patr*0476 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0477         use OAD_cp
                0478         implicit none
                0479         double precision, dimension(:,:,:,:) :: x
8702af1f36 Patr*0480 C $OpenAD$ END DECLS
9c5e597b65 Jean*0481         read(unit=cp_io_unit) x
8702af1f36 Patr*0482 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0483         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1)
8702af1f36 Patr*0484 #endif
                0485       end subroutine 
                0486 
9c5e597b65 Jean*0487       subroutine cp_arg_store_real_four_tensor_a(x)
8702af1f36 Patr*0488 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0489         use OAD_active
                0490         use OAD_cp
                0491         implicit none
                0492         type(active), dimension(:,:,:,:) :: x
8702af1f36 Patr*0493 C $OpenAD$ END DECLS
                0494 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0495         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1)%v
8702af1f36 Patr*0496 #endif
9c5e597b65 Jean*0497         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0498       end subroutine 
                0499 
                0500       subroutine cp_arg_restore_real_four_tensor_a(x)
                0501 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0502         use OAD_active
                0503         use OAD_cp
                0504         implicit none
                0505         type(active), dimension(:,:,:,:) :: x
8702af1f36 Patr*0506 C $OpenAD$ END DECLS
9c5e597b65 Jean*0507         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0508 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0509         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1)%v
7d8fa04723 Jean*0510 #endif 
                0511       end subroutine
8702af1f36 Patr*0512 
                0513       subroutine cp_arg_store_real_five_tensor(x)
                0514 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0515         use OAD_cp
                0516         implicit none
                0517         double precision, dimension(:,:,:,:,:) :: x
8702af1f36 Patr*0518 C $OpenAD$ END DECLS
                0519 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0520         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1,1)
8702af1f36 Patr*0521 #endif
9c5e597b65 Jean*0522         write(unit=cp_io_unit) x
8702af1f36 Patr*0523       end subroutine 
                0524 
9c5e597b65 Jean*0525       subroutine cp_arg_restore_real_five_tensor(x)
8702af1f36 Patr*0526 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0527         use OAD_cp
                0528         implicit none
                0529         double precision, dimension(:,:,:,:,:) :: x
8702af1f36 Patr*0530 C $OpenAD$ END DECLS
9c5e597b65 Jean*0531         read(unit=cp_io_unit) x
8702af1f36 Patr*0532 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0533         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1,1)
8702af1f36 Patr*0534 #endif
                0535       end subroutine 
                0536 
9c5e597b65 Jean*0537       subroutine cp_arg_store_real_five_tensor_a(x)
8702af1f36 Patr*0538 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0539         use OAD_active
                0540         use OAD_cp
                0541         implicit none
                0542         type(active), dimension(:,:,:,:,:) :: x
8702af1f36 Patr*0543 C $OpenAD$ END DECLS
                0544 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0545         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1,1)%v
8702af1f36 Patr*0546 #endif
9c5e597b65 Jean*0547         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0548       end subroutine 
                0549 
                0550       subroutine cp_arg_restore_real_five_tensor_a(x)
                0551 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0552         use OAD_active
                0553         use OAD_cp
                0554         implicit none
                0555         type(active), dimension(:,:,:,:,:) :: x
8702af1f36 Patr*0556 C $OpenAD$ END DECLS
9c5e597b65 Jean*0557         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0558 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0559         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1,1)%v
8702af1f36 Patr*0560 #endif
                0561       end subroutine 
                0562 
7d8fa04723 Jean*0563       subroutine cp_arg_store_real_six_tensor(x)
8702af1f36 Patr*0564 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0565         use OAD_cp
                0566         implicit none
                0567         double precision, dimension(:,:,:,:,:,:) :: x
8702af1f36 Patr*0568 C $OpenAD$ END DECLS
                0569 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0570         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1,1,1)
8702af1f36 Patr*0571 #endif
9c5e597b65 Jean*0572         write(unit=cp_io_unit) x
8702af1f36 Patr*0573       end subroutine 
                0574 
7d8fa04723 Jean*0575       subroutine cp_arg_restore_real_six_tensor(x)
8702af1f36 Patr*0576 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0577         use OAD_cp
                0578         implicit none
                0579         double precision, dimension(:,:,:,:,:,:) :: x
8702af1f36 Patr*0580 C $OpenAD$ END DECLS
9c5e597b65 Jean*0581         read(unit=cp_io_unit) x
8702af1f36 Patr*0582 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0583         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1,1,1)
8702af1f36 Patr*0584 #endif
                0585       end subroutine 
                0586 
7d8fa04723 Jean*0587       subroutine cp_arg_store_real_six_tensor_a(x)
8702af1f36 Patr*0588 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0589         use OAD_active
                0590         use OAD_cp
                0591         implicit none
                0592         type(active), dimension(:,:,:,:,:,:) :: x
8702af1f36 Patr*0593 C $OpenAD$ END DECLS
                0594 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0595         !write(standardmessageunit,*)'OAD: cp write x ', x(1,1,1,1,1,1)%v
8702af1f36 Patr*0596 #endif
7d8fa04723 Jean*0597         write(unit=cp_io_unit) x%v
8702af1f36 Patr*0598       end subroutine 
                0599 
7d8fa04723 Jean*0600       subroutine cp_arg_restore_real_six_tensor_a(x)
8702af1f36 Patr*0601 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0602         use OAD_active
                0603         use OAD_cp
                0604         implicit none
                0605         type(active), dimension(:,:,:,:,:,:) :: x
8702af1f36 Patr*0606 C $OpenAD$ END DECLS
7d8fa04723 Jean*0607         read(unit=cp_io_unit) x%v
8702af1f36 Patr*0608 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0609         !write(standardmessageunit,*)'OAD: cp read x ', x(1,1,1,1,1,1)%v
8702af1f36 Patr*0610 #endif
                0611       end subroutine 
                0612 
7d8fa04723 Jean*0613 C integers -----------------------------------------------------
                0614       subroutine cp_arg_store_integer_scalar(i)
8702af1f36 Patr*0615 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0616         use OAD_cp
                0617         implicit none
                0618         integer :: i
8702af1f36 Patr*0619 C $OpenAD$ END DECLS
                0620 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0621         !write(standardmessageunit,*)'OAD: cp write i ', i
8702af1f36 Patr*0622 #endif
7d8fa04723 Jean*0623         write(unit=cp_io_unit) i
8702af1f36 Patr*0624       end subroutine 
                0625 
7d8fa04723 Jean*0626       subroutine cp_arg_restore_integer_scalar(i)
8702af1f36 Patr*0627 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0628         use OAD_cp
                0629         implicit none
                0630         integer :: i
8702af1f36 Patr*0631 C $OpenAD$ END DECLS
7d8fa04723 Jean*0632         read(unit=cp_io_unit) i
8702af1f36 Patr*0633 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0634         !write(standardmessageunit,*)'OAD: cp read i ', i
8702af1f36 Patr*0635 #endif
                0636       end subroutine 
                0637 
7d8fa04723 Jean*0638       subroutine cp_arg_store_integer_vector(i)
8702af1f36 Patr*0639 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0640         use OAD_cp
                0641         implicit none
                0642         integer, dimension(:) :: i
8702af1f36 Patr*0643 C $OpenAD$ END DECLS
                0644 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0645         !write(standardmessageunit,*)'OAD: cp write i ', i(1)
8702af1f36 Patr*0646 #endif
7d8fa04723 Jean*0647         write(unit=cp_io_unit) i
8702af1f36 Patr*0648       end subroutine 
                0649 
7d8fa04723 Jean*0650       subroutine cp_arg_restore_integer_vector(i)
8702af1f36 Patr*0651 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0652         use OAD_cp
                0653         implicit none
                0654         integer, dimension(:) :: i
8702af1f36 Patr*0655 C $OpenAD$ END DECLS
7d8fa04723 Jean*0656         read(unit=cp_io_unit) i
8702af1f36 Patr*0657 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0658         !write(standardmessageunit,*)'OAD: cp read i ', i(1)
8702af1f36 Patr*0659 #endif
                0660       end subroutine 
                0661 
7d8fa04723 Jean*0662       subroutine cp_arg_store_integer_matrix(i)
8702af1f36 Patr*0663 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0664         use OAD_cp
                0665         implicit none
                0666         integer, dimension(:,:) :: i
8702af1f36 Patr*0667 C $OpenAD$ END DECLS
                0668 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0669         !write(standardmessageunit,*)'OAD: cp write i ', i(1,1)
8702af1f36 Patr*0670 #endif
7d8fa04723 Jean*0671         write(unit=cp_io_unit) i
8702af1f36 Patr*0672       end subroutine 
                0673 
7d8fa04723 Jean*0674       subroutine cp_arg_restore_integer_matrix(i)
8702af1f36 Patr*0675 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0676         use OAD_cp
                0677         implicit none
                0678         integer, dimension(:,:) :: i
8702af1f36 Patr*0679 C $OpenAD$ END DECLS
7d8fa04723 Jean*0680         read(unit=cp_io_unit) i
8702af1f36 Patr*0681 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0682         !write(standardmessageunit,*)'OAD: cp read i ', i(1,1)
8702af1f36 Patr*0683 #endif
                0684       end subroutine 
                0685 
7d8fa04723 Jean*0686       subroutine cp_arg_store_integer_three_tensor(i)
8702af1f36 Patr*0687 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0688         use OAD_cp
                0689         implicit none
                0690         integer, dimension(:,:,:) :: i
8702af1f36 Patr*0691 C $OpenAD$ END DECLS
                0692 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0693         !write(standardmessageunit,*)'OAD: cp write i ', i(1,1,1)
8702af1f36 Patr*0694 #endif
7d8fa04723 Jean*0695         write(unit=cp_io_unit) i
8702af1f36 Patr*0696       end subroutine 
                0697 
7d8fa04723 Jean*0698       subroutine cp_arg_restore_integer_three_tensor(i)
8702af1f36 Patr*0699 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0700         use OAD_cp
                0701         implicit none
                0702         integer, dimension(:,:,:) :: i
8702af1f36 Patr*0703 C $OpenAD$ END DECLS
7d8fa04723 Jean*0704         read(unit=cp_io_unit) i
8702af1f36 Patr*0705 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0706         !write(standardmessageunit,*)'OAD: cp read i ', i(1,1,1)
7d8fa04723 Jean*0707 #endif
                0708       end subroutine 
                0709 
                0710       subroutine cp_arg_store_integer_four_tensor(i)
                0711 C $OpenAD$ INLINE DECLS
                0712         use OAD_cp
                0713         implicit none
                0714         integer, dimension(:,:,:,:) :: i
                0715 C $OpenAD$ END DECLS
                0716 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0717         !write(standardmessageunit,*)'OAD: cp write i ', i(1,1,1,1)
7d8fa04723 Jean*0718 #endif
                0719         write(unit=cp_io_unit) i
                0720       end subroutine 
                0721 
                0722       subroutine cp_arg_restore_integer_four_tensor(i)
                0723 C $OpenAD$ INLINE DECLS
                0724         use OAD_cp
                0725         implicit none
                0726         integer, dimension(:,:,:,:) :: i
                0727 C $OpenAD$ END DECLS
                0728         read(unit=cp_io_unit) i
                0729 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0730         !write(standardmessageunit,*)'OAD: cp read i ', i(1,1,1,1)
7d8fa04723 Jean*0731 #endif
                0732       end subroutine 
                0733 
                0734       subroutine cp_arg_store_integer_five_tensor(i)
                0735 C $OpenAD$ INLINE DECLS
                0736         use OAD_cp
                0737         implicit none
                0738         integer, dimension(:,:,:,:,:) :: i
                0739 C $OpenAD$ END DECLS
                0740 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0741         !write(standardmessageunit,*)'OAD: cp write i ', i(1,1,1,1,1)
7d8fa04723 Jean*0742 #endif
                0743         write(unit=cp_io_unit) i
                0744       end subroutine 
                0745 
                0746       subroutine cp_arg_restore_integer_five_tensor(i)
                0747 C $OpenAD$ INLINE DECLS
                0748         use OAD_cp
                0749         implicit none
                0750         integer, dimension(:,:,:,:,:) :: i
                0751 C $OpenAD$ END DECLS
                0752         read (unit=cp_io_unit) i
                0753 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0754         !write(standardmessageunit,*)'OAD: cp read i ', i(1,1,1,1,1)
8702af1f36 Patr*0755 #endif
                0756       end subroutine 
                0757 
                0758 C strings  -----------------------------------------------------
7d8fa04723 Jean*0759       subroutine cp_arg_store_string_scalar(s)
8702af1f36 Patr*0760 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0761       use OAD_cp
8702af1f36 Patr*0762       implicit none
7d8fa04723 Jean*0763       character*(80) :: s
8702af1f36 Patr*0764 C $OpenAD$ END DECLS 
                0765 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0766         !write(standardmessageunit,*)'OAD: cp write s ', s
8702af1f36 Patr*0767 #endif
7d8fa04723 Jean*0768         write(unit=cp_io_unit) s
8702af1f36 Patr*0769       end subroutine 
                0770       
7d8fa04723 Jean*0771       subroutine cp_arg_restore_string_scalar(s)
8702af1f36 Patr*0772 C $OpenAD$ INLINE DECLS
9c5e597b65 Jean*0773       use OAD_cp
8702af1f36 Patr*0774       implicit none
7d8fa04723 Jean*0775       character*(80) :: s
8702af1f36 Patr*0776 C $OpenAD$ END DECLS
7d8fa04723 Jean*0777         read (unit=cp_io_unit) s
8702af1f36 Patr*0778 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0779         !write(standardmessageunit,*)'OAD: cp read s ', s
8702af1f36 Patr*0780 #endif
                0781       end subroutine 
                0782 
                0783 C bools  -----------------------------------------------------
7d8fa04723 Jean*0784       subroutine cp_arg_store_bool_scalar(b)
8702af1f36 Patr*0785 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0786         use OAD_cp
                0787         implicit none
                0788         logical :: b
8702af1f36 Patr*0789 C $OpenAD$ END DECLS
                0790 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0791         !write(standardmessageunit,*)'OAD: cp write b ', b
8702af1f36 Patr*0792 #endif
7d8fa04723 Jean*0793         write(unit=cp_io_unit) b
8702af1f36 Patr*0794       end subroutine 
                0795 
7d8fa04723 Jean*0796       subroutine cp_arg_restore_bool_scalar(b)
8702af1f36 Patr*0797 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0798         use OAD_cp
                0799         implicit none
                0800         logical :: b
8702af1f36 Patr*0801 C $OpenAD$ END DECLS
7d8fa04723 Jean*0802         read (unit=cp_io_unit) b
8702af1f36 Patr*0803 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0804         !write(standardmessageunit,*)'OAD: cp read b ', b
8702af1f36 Patr*0805 #endif
                0806       end subroutine 
9c5e597b65 Jean*0807 
7d8fa04723 Jean*0808       subroutine cp_arg_store_bool_matrix(b)
9c5e597b65 Jean*0809 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0810         use OAD_cp
                0811         implicit none
                0812         logical, dimension(:,:) :: b
9c5e597b65 Jean*0813 C $OpenAD$ END DECLS
766df7f670 Oliv*0814 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0815         !write(standardmessageunit,*)'OAD: cp write b ', b(1,1)
9c5e597b65 Jean*0816 #endif
7d8fa04723 Jean*0817         write(unit=cp_io_unit) b
9c5e597b65 Jean*0818       end subroutine 
                0819 
7d8fa04723 Jean*0820       subroutine cp_arg_restore_bool_matrix(b)
9c5e597b65 Jean*0821 C $OpenAD$ INLINE DECLS
7d8fa04723 Jean*0822         use OAD_cp
                0823         implicit none
                0824         logical, dimension(:,:) :: b
9c5e597b65 Jean*0825 C $OpenAD$ END DECLS
7d8fa04723 Jean*0826         read(unit=cp_io_unit) b
9c5e597b65 Jean*0827 #ifdef OAD_DEBUG_CP
279dc77b07 Patr*0828         !write(standardmessageunit,*)'OAD: cp read b ', b(1,1)
9c5e597b65 Jean*0829 #endif
279dc77b07 Patr*0830       end subroutine
                0831 
                0832 C adjoints of active reals ----------------------------------
                0833       subroutine cp_arg_store_real_scalar_a_d(x)
                0834 C $OpenAD$ INLINE DECLS
                0835         use OAD_active
                0836         use OAD_cp
                0837         implicit none
                0838         type(active) :: x
                0839 C $OpenAD$ END DECLS
                0840 #ifdef OAD_DEBUG_CP
                0841         !write(standardmessageunit,*)'OAD: cp write x%d ', x%d
                0842 #endif
                0843         write(unit=cp_io_unit) x%d
                0844       end subroutine 
                0845 
                0846       subroutine cp_arg_restore_real_scalar_a_d(x)
                0847 C $OpenAD$ INLINE DECLS
                0848         use OAD_active
                0849         use OAD_cp
                0850         implicit none
                0851         type(active) :: x
                0852 C $OpenAD$ END DECLS
                0853         read(unit=cp_io_unit) x%d
                0854 #ifdef OAD_DEBUG_CP
                0855         !write(standardmessageunit,*)'OAD: cp read x%d ', x%d
                0856 #endif
                0857       end subroutine 
                0858 
                0859       subroutine cp_arg_store_real_vector_a_d(x)
                0860 C $OpenAD$ INLINE DECLS
                0861         use OAD_active
                0862         use OAD_cp
                0863         implicit none
                0864         type(active), dimension(:) :: x
                0865 C $OpenAD$ END DECLS
                0866 #ifdef OAD_DEBUG_CP
                0867         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1)%d
                0868 #endif
                0869         write(unit=cp_io_unit) x%d
                0870       end subroutine 
                0871 
                0872       subroutine cp_arg_restore_real_vector_a_d(x)
                0873 C $OpenAD$ INLINE DECLS
                0874         use OAD_active
                0875         use OAD_cp
                0876         implicit none
                0877         type(active), dimension(:) :: x
                0878 C $OpenAD$ END DECLS
                0879         read(unit=cp_io_unit) x%d
                0880 #ifdef OAD_DEBUG_CP
                0881         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1)%d
                0882 #endif
                0883       end subroutine  
                0884 
                0885       subroutine cp_arg_store_real_matrix_a_d(x)
                0886 C $OpenAD$ INLINE DECLS
                0887         use OAD_active
                0888         use OAD_cp
                0889         implicit none
                0890         type(active), dimension(:,:) :: x
                0891 C $OpenAD$ END DECLS
                0892 #ifdef OAD_DEBUG_CP
                0893         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1,1)%d
                0894 #endif
                0895         write(unit=cp_io_unit) x%d
                0896       end subroutine 
                0897 
                0898       subroutine cp_arg_restore_real_matrix_a_d(x)
                0899 C $OpenAD$ INLINE DECLS
                0900         use OAD_active
                0901         use OAD_cp
                0902         implicit none
                0903         type(active), dimension(:,:) :: x
                0904 C $OpenAD$ END DECLS
                0905         read(unit=cp_io_unit) x%d
                0906 #ifdef OAD_DEBUG_CP
                0907         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1,1)%d
                0908 #endif
                0909       end subroutine 
                0910 
                0911       subroutine cp_arg_store_real_three_tensor_a_d(x)
                0912 C $OpenAD$ INLINE DECLS
                0913         use OAD_active
                0914         use OAD_cp
                0915         implicit none
                0916         type(active), dimension(:,:,:) :: x
                0917 C $OpenAD$ END DECLS
                0918 #ifdef OAD_DEBUG_CP
                0919         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1,1,1)%d
                0920 #endif
                0921         write(unit=cp_io_unit) x%d
9c5e597b65 Jean*0922       end subroutine 
279dc77b07 Patr*0923 
                0924       subroutine cp_arg_restore_real_three_tensor_a_d(x)
                0925 C $OpenAD$ INLINE DECLS
                0926         use OAD_active
                0927         use OAD_cp
                0928         implicit none
                0929         type(active), dimension(:,:,:) :: x
                0930 C $OpenAD$ END DECLS
                0931         read(unit=cp_io_unit) x%d
                0932 #ifdef OAD_DEBUG_CP
                0933         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1,1,1)%d
                0934 #endif
                0935       end subroutine 
                0936 
                0937       subroutine cp_arg_store_real_four_tensor_a_d(x)
                0938 C $OpenAD$ INLINE DECLS
                0939         use OAD_active
                0940         use OAD_cp
                0941         implicit none
                0942         type(active), dimension(:,:,:,:) :: x
                0943 C $OpenAD$ END DECLS
                0944 #ifdef OAD_DEBUG_CP
                0945         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1,1,1,1)%d
                0946 #endif
                0947         write(unit=cp_io_unit) x%d
                0948       end subroutine 
                0949 
                0950       subroutine cp_arg_restore_real_four_tensor_a_d(x)
                0951 C $OpenAD$ INLINE DECLS
                0952         use OAD_active
                0953         use OAD_cp
                0954         implicit none
                0955         type(active), dimension(:,:,:,:) :: x
                0956 C $OpenAD$ END DECLS
                0957         read(unit=cp_io_unit) x%d
                0958 #ifdef OAD_DEBUG_CP
                0959         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1,1,1,1)%d
                0960 #endif 
                0961       end subroutine
                0962 
                0963       subroutine cp_arg_store_real_five_tensor_a_d(x)
                0964 C $OpenAD$ INLINE DECLS
                0965         use OAD_active
                0966         use OAD_cp
                0967         implicit none
                0968         type(active), dimension(:,:,:,:,:) :: x
                0969 C $OpenAD$ END DECLS
                0970 #ifdef OAD_DEBUG_CP
                0971         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1,1,1,1,1)%d
                0972 #endif
                0973         write(unit=cp_io_unit) x%d
                0974       end subroutine 
                0975 
                0976       subroutine cp_arg_restore_real_five_tensor_a_d(x)
                0977 C $OpenAD$ INLINE DECLS
                0978         use OAD_active
                0979         use OAD_cp
                0980         implicit none
                0981         type(active), dimension(:,:,:,:,:) :: x
                0982 C $OpenAD$ END DECLS
                0983         read(unit=cp_io_unit) x%d
                0984 #ifdef OAD_DEBUG_CP
                0985         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1,1,1,1,1)%d
                0986 #endif
                0987       end subroutine 
                0988 
                0989       subroutine cp_arg_store_real_six_tensor_a_d(x)
                0990 C $OpenAD$ INLINE DECLS
                0991         use OAD_active
                0992         use OAD_cp
                0993         implicit none
                0994         type(active), dimension(:,:,:,:,:,:) :: x
                0995 C $OpenAD$ END DECLS
                0996 #ifdef OAD_DEBUG_CP
                0997         !write(standardmessageunit,*)'OAD: cp write x%d ', x(1,1,1,1,1,1)%d
                0998 #endif
                0999         write(unit=cp_io_unit) x%d
                1000       end subroutine 
                1001 
                1002       subroutine cp_arg_restore_real_six_tensor_a_d(x)
                1003 C $OpenAD$ INLINE DECLS
                1004         use OAD_active
                1005         use OAD_cp
                1006         implicit none
                1007         type(active), dimension(:,:,:,:,:,:) :: x
                1008 C $OpenAD$ END DECLS
                1009         read(unit=cp_io_unit) x%d
                1010 #ifdef OAD_DEBUG_CP
                1011         !write(standardmessageunit,*)'OAD: cp read x%d ', x(1,1,1,1,1,1)%d
                1012 #endif
                1013       end subroutine