#include "cctk.h" #include "cctk_Functions.h" #include "cctk_Parameters.h" subroutine deriv_gf_8_4_opt ( var, ni, nj, nk, dir, bb, gsize, delta, dvar ) implicit none DECLARE_CCTK_FUNCTIONS DECLARE_CCTK_PARAMETERS CCTK_REAL, parameter :: zero = 0.0 integer, parameter :: wp = kind(zero) CCTK_INT, intent(IN) :: ni, nj, nk CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: var CCTK_INT, intent(IN) :: dir CCTK_INT, intent(IN) :: bb(2) CCTK_INT, intent(IN) :: gsize CCTK_REAL, intent(IN) :: delta CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar CCTK_REAL, dimension(4), save :: a CCTK_REAL, dimension(12,8), save :: q CCTK_REAL :: idel CCTK_INT :: il, ir, jl, jr, kl, kr, j, k logical, save :: first = .true. if ( first ) then a(1) = 4.0_wp/5.0_wp; a(2) = -1.0_wp/5.0_wp a(3) = 4.0_wp/105.0_wp; a(4) = -1.0_wp/280.0_wp q(1,1) = -1.695543604431898508749855654248370812054_wp q(2,1) = 2.244525109969933844775276252674718713114_wp q(3,1) = -0.002163206607021876385102025416828692597556_wp q(4,1) = -0.8963677947307838167920190976769921451579_wp q(5,1) = 0.2272938673777916084251861567415263328923_wp q(6,1) = 0.1564962758424078514454013191006671546757_wp q(7,1) = 0.003067629343740613524557473056482141287158_wp q(8,1) = -0.03730827676416971624344442423120269215938_wp q(9,1) = 0.0_wp q(10,1) = 0.0_wp q(11,1) = 0.0_wp q(12,1) = 0.0_wp q(1,2) = -0.4338209217401506177055540526837828066711_wp q(2,2) = 0.0_wp q(3,2) = 0.06631961192627293470305188523936353462439_wp q(4,2) = 0.6021342012739069461791324177355805603656_wp q(5,2) = -0.2213339967864188453859900127368926875520_wp q(6,2) = -0.001759980003751516130708214237001705688807_wp q(7,2) = -0.03586646217917856587352426166797879385310_wp q(8,2) = 0.02432754750931966421359223835071189877503_wp q(9,2) = 0.0_wp q(10,2) = 0.0_wp q(11,2) = 0.0_wp q(12,2) = 0.0_wp q(1,3) = 0.002477771724790106958560398469919805908522_wp q(2,3) = -0.3930241559935857537756812928554282732813_wp q(3,3) = 0.0_wp q(4,3) = -0.6080263371305341035810612202824357742461_wp q(5,3) = 2.372533124692864412954923458457782022478_wp q(6,3) = -2.126112217612194556872952028715623246817_wp q(7,3) = 0.9075309435559160631495435310307745215131_wp q(8,3) = -0.1553791292372561688333328461049890555549_wp q(9,3) = 0.0_wp q(10,3) = 0.0_wp q(11,3) = 0.0_wp q(12,3) = 0.0_wp q(1,4) = 0.1470043328582935680963279007411074464266_wp q(2,4) = -0.5109179516689331400658365102850902594996_wp q(3,4) = 0.08705685833207777685654755900085198684605_wp q(4,4) = 0.0_wp q(5,4) = -0.1306945333157231839535552435788086367322_wp q(6,4) = 0.5981933016362239542598518794410035592720_wp q(7,4) = -0.2031099149801233772210520666270517105660_wp q(8,4) = 0.01246790713818440202771648130798761425319_wp q(9,4) = 0.0_wp q(10,4) = 0.0_wp q(11,4) = 0.0_wp q(12,4) = 0.0_wp q(1,5) = -0.1624073990846984662267508265053107632238_wp q(2,5) = 0.8182390368133059538132603839842499379861_wp q(3,5) = -1.480018301574606037840376638130713134769_wp q(4,5) = 0.5694185675497882973361521309100738568707_wp q(5,5) = 0.0_wp q(6,5) = 0.02929508293831017534661096502884570435894_wp q(7,5) = 0.2825910005984090924748504462528740123097_wp q(8,5) = -0.04846434332860791750685068920563907354793_wp q(9,5) = -0.008653643911901097396895772334380539984709_wp q(10,5) = 0.0_wp q(11,5) = 0.0_wp q(12,5) = 0.0_wp q(1,6) = -0.03609686950604370828405582087760384548297_wp q(2,6) = 0.002100328577309696555612805397001149641337_wp q(3,6) = 0.4281425817419204360479874008765051769418_wp q(4,6) = -0.8413238238875046736839550701643761633826_wp q(5,6) = -0.009456755727630000971085851751473603919416_wp q(6,6) = 0.0_wp q(7,6) = 0.5172389789041733771234394082779034212042_wp q(8,6) = -0.08760813565107717873901807278223424202774_wp q(9,6) = 0.02979718129528502284256573906127239395909_wp q(10,6) = -0.002793485746432970891490538036994286933664_wp q(11,6) = 0.0_wp q(12,6) = 0.0_wp q(1,7) = -0.0009797678134978722516935350416937185004514_wp q(2,7) = 0.05926834509975463070197111976550139351983_wp q(3,7) = -0.2530570463899371286637621744353665227114_wp q(4,7) = 0.3955555826583941989223787901885502942688_wp q(5,7) = -0.1263166266018192756531792392597193430222_wp q(6,7) = -0.7162192643577544899896537844517283192029_wp q(7,7) = 0.0_wp q(8,7) = 0.8209721963136350137518635528607308559584_wp q(9,7) = -0.2166153552278720352907291696202456084323_wp q(10,7) = 0.04126006766245181624585317516576106827282_wp q(11,7) = -0.003868131343354857773048735171790100150577_wp q(12,7) = 0.0_wp q(1,8) = 0.01090012273307913185972171872891927027272_wp q(2,8) = -0.03677379943661633440187210478144113487019_wp q(3,8) = 0.03963287609450569641558898819403594845401_wp q(4,8) = -0.02221139656912423686782339404035534965598_wp q(5,8) = 0.01981665906734246925389947980969432794418_wp q(6,8) = 0.1109698768905366566973791789571759073485_wp q(7,8) = -0.7509903604688148129224205834189298636300_wp q(8,8) = 0.0_wp q(9,8) = 0.7926019635554773751160111698352821779736_wp q(10,8) = -0.1981504908888693437790027924588205444934_wp q(11,8) = 0.03774295064549892262457196046834677037969_wp q(12,8) = -0.003538401623015523996053621293907509723096_wp first = .false. end if idel = 1.0_wp / delta direction: select case (dir) case (0) direction if ( bb(1) == 0 ) then il = 1 + gsize else dvar(1,:,:) = ( q(1,1) * var(1,:,:) + q(2,1) * var(2,:,:) + & q(3,1) * var(3,:,:) + q(4,1) * var(4,:,:) + & q(5,1) * var(5,:,:) + q(6,1) * var(6,:,:) + & q(7,1) * var(7,:,:) + q(8,1) * var(8,:,:) ) * idel dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(3,2) * var(3,:,:) + & q(4,2) * var(4,:,:) + q(5,2) * var(5,:,:) + & q(6,2) * var(6,:,:) + q(7,2) * var(7,:,:) + & q(8,2) * var(8,:,:) ) * idel dvar(3,:,:) = ( q(1,3) * var(1,:,:) + q(2,3) * var(2,:,:) + & q(4,3) * var(4,:,:) + q(5,3) * var(5,:,:) + & q(6,3) * var(6,:,:) + q(7,3) * var(7,:,:) + & q(8,3) * var(8,:,:) ) * idel dvar(4,:,:) = ( q(1,4) * var(1,:,:) + q(2,4) * var(2,:,:) + & q(3,4) * var(3,:,:) + q(5,4) * var(5,:,:) + & q(6,4) * var(6,:,:) + q(7,4) * var(7,:,:) + & q(8,4) * var(8,:,:) ) * idel dvar(5,:,:) = ( q(1,5) * var(1,:,:) + q(2,5) * var(2,:,:) + & q(3,5) * var(3,:,:) + q(4,5) * var(4,:,:) + & q(6,5) * var(6,:,:) + q(7,5) * var(7,:,:) + & q(8,5) * var(8,:,:) + q(9,5) * var(9,:,:) ) * idel dvar(6,:,:) = ( q(1,6) * var(1,:,:) + q(2,6) * var(2,:,:) + & q(3,6) * var(3,:,:) + q(4,6) * var(4,:,:) + & q(5,6) * var(5,:,:) + q(7,6) * var(7,:,:) + & q(8,6) * var(8,:,:) + q(9,6) * var(9,:,:) + & q(10,6) * var(10,:,:) ) * idel dvar(7,:,:) = ( q(1,7) * var(1,:,:) + q(2,7) * var(2,:,:) + & q(3,7) * var(3,:,:) + q(4,7) * var(4,:,:) + & q(5,7) * var(5,:,:) + q(6,7) * var(6,:,:) + & q(8,7) * var(8,:,:) + q(9,7) * var(9,:,:) + & q(10,7) * var(10,:,:) + q(11,7) * var(11,:,:) ) * idel dvar(8,:,:) = ( q(1,8) * var(1,:,:) + q(2,8) * var(2,:,:) + & q(3,8) * var(3,:,:) + q(4,8) * var(4,:,:) + & q(5,8) * var(5,:,:) + q(6,8) * var(6,:,:) + & q(7,8) * var(7,:,:) + q(9,8) * var(9,:,:) + & q(10,8) * var(10,:,:) + q(11,8) * var(11,:,:) + & q(12,8) * var(12,:,:) ) * idel il = 9 end if if ( bb(2) == 0 ) then ir = ni - gsize else dvar(ni,:,:) = - ( q(1,1) * var(ni,:,:) + q(2,1) * var(ni-1,:,:) + & q(3,1) * var(ni-2,:,:) + q(4,1) * var(ni-3,:,:) + & q(5,1) * var(ni-4,:,:) + q(6,1) * var(ni-5,:,:) + & q(7,1) * var(ni-6,:,:) + & q(8,1) * var(ni-7,:,:) ) * idel dvar(ni-1,:,:) = - ( q(1,2) * var(ni,:,:) + q(3,2) * var(ni-2,:,:) + & q(4,2) * var(ni-3,:,:) + q(5,2) * var(ni-4,:,:) + & q(6,2) * var(ni-5,:,:) + q(7,2) * var(ni-6,:,:) + & q(8,2) * var(ni-7,:,:) ) * idel dvar(ni-2,:,:) = - ( q(1,3) * var(ni,:,:) + q(2,3) * var(ni-1,:,:) + & q(4,3) * var(ni-3,:,:) + q(5,3) * var(ni-4,:,:) + & q(6,3) * var(ni-5,:,:) + q(7,3) * var(ni-6,:,:) + & q(8,3) * var(ni-7,:,:) ) * idel dvar(ni-3,:,:) = - ( q(1,4) * var(ni,:,:) + q(2,4) * var(ni-1,:,:) + & q(3,4) * var(ni-2,:,:) + q(5,4) * var(ni-4,:,:) + & q(6,4) * var(ni-5,:,:) + q(7,4) * var(ni-6,:,:) + & q(8,4) * var(ni-7,:,:) ) * idel dvar(ni-4,:,:) = - ( q(1,5) * var(ni,:,:) + q(2,5) * var(ni-1,:,:) + & q(3,5) * var(ni-2,:,:) + q(4,5) * var(ni-3,:,:) + & q(6,5) * var(ni-5,:,:) + q(7,5) * var(ni-6,:,:) + & q(8,5) * var(ni-7,:,:) + & q(9,5) * var(ni-8,:,:) ) * idel dvar(ni-5,:,:) = - ( q(1,6) * var(ni,:,:) + q(2,6) * var(ni-1,:,:) + & q(3,6) * var(ni-2,:,:) + q(4,6) * var(ni-3,:,:) + & q(5,6) * var(ni-4,:,:) + q(7,6) * var(ni-6,:,:) + & q(8,6) * var(ni-7,:,:) + q(9,6) * var(ni-8,:,:) + & q(10,6) * var(ni-9,:,:) ) * idel dvar(ni-6,:,:) = - ( q(1,7) * var(ni,:,:) + q(2,7) * var(ni-1,:,:) + & q(3,7) * var(ni-2,:,:) + q(4,7) * var(ni-3,:,:) + & q(5,7) * var(ni-4,:,:) + q(6,7) * var(ni-5,:,:) + & q(8,7) * var(ni-7,:,:) + q(9,7) * var(ni-8,:,:) + & q(10,7) * var(ni-9,:,:) + & q(11,7) * var(ni-10,:,:) ) * idel dvar(ni-7,:,:) = - ( q(1,8) * var(ni,:,:) + q(2,8) * var(ni-1,:,:) + & q(3,8) * var(ni-2,:,:) + q(4,8) * var(ni-3,:,:) + & q(5,8) * var(ni-4,:,:) + q(6,8) * var(ni-5,:,:) + & q(7,8) * var(ni-6,:,:) + q(9,8) * var(ni-8,:,:) + & q(10,8) * var(ni-9,:,:) + & q(11,8) * var(ni-10,:,:) + & q(12,8) * var(ni-11,:,:) ) * idel ir = ni - 8 end if dvar(il:ir,:,:) = ( a(1) * ( var(il+1:ir+1,:,:) - & var(il-1:ir-1,:,:) ) + & a(2) * ( var(il+2:ir+2,:,:) - & var(il-2:ir-2,:,:) ) + & a(3) * ( var(il+3:ir+3,:,:) - & var(il-3:ir-3,:,:) ) + & a(4) * ( var(il+4:ir+4,:,:) - & var(il-4:ir-4,:,:) ) ) * idel case (1) direction if ( zero_derivs_y /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then jl = 1 + gsize else dvar(:,1,:) = ( q(1,1) * var(:,1,:) + q(2,1) * var(:,2,:) + & q(3,1) * var(:,3,:) + q(4,1) * var(:,4,:) + & q(5,1) * var(:,5,:) + q(6,1) * var(:,6,:) + & q(7,1) * var(:,7,:) + q(8,1) * var(:,8,:) ) * idel dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(3,2) * var(:,3,:) + & q(4,2) * var(:,4,:) + q(5,2) * var(:,5,:) + & q(6,2) * var(:,6,:) + q(7,2) * var(:,7,:) + & q(8,2) * var(:,8,:) ) * idel dvar(:,3,:) = ( q(1,3) * var(:,1,:) + q(2,3) * var(:,2,:) + & q(4,3) * var(:,4,:) + q(5,3) * var(:,5,:) + & q(6,3) * var(:,6,:) + q(7,3) * var(:,7,:) + & q(8,3) * var(:,8,:) ) * idel dvar(:,4,:) = ( q(1,4) * var(:,1,:) + q(2,4) * var(:,2,:) + & q(3,4) * var(:,3,:) + q(5,4) * var(:,5,:) + & q(6,4) * var(:,6,:) + q(7,4) * var(:,7,:) + & q(8,4) * var(:,8,:) ) * idel dvar(:,5,:) = ( q(1,5) * var(:,1,:) + q(2,5) * var(:,2,:) + & q(3,5) * var(:,3,:) + q(4,5) * var(:,4,:) + & q(6,5) * var(:,6,:) + q(7,5) * var(:,7,:) + & q(8,5) * var(:,8,:) + q(9,5) * var(:,9,:) ) * idel dvar(:,6,:) = ( q(1,6) * var(:,1,:) + q(2,6) * var(:,2,:) + & q(3,6) * var(:,3,:) + q(4,6) * var(:,4,:) + & q(5,6) * var(:,5,:) + q(7,6) * var(:,7,:) + & q(8,6) * var(:,8,:) + q(9,6) * var(:,9,:) + & q(10,6) * var(:,10,:) ) * idel dvar(:,7,:) = ( q(1,7) * var(:,1,:) + q(2,7) * var(:,2,:) + & q(3,7) * var(:,3,:) + q(4,7) * var(:,4,:) + & q(5,7) * var(:,5,:) + q(6,7) * var(:,6,:) + & q(8,7) * var(:,8,:) + q(9,7) * var(:,9,:) + & q(10,7) * var(:,10,:) + q(11,7) * var(:,11,:) ) * idel dvar(:,8,:) = ( q(1,8) * var(:,1,:) + q(2,8) * var(:,2,:) + & q(3,8) * var(:,3,:) + q(4,8) * var(:,4,:) + & q(5,8) * var(:,5,:) + q(6,8) * var(:,6,:) + & q(7,8) * var(:,7,:) + q(9,8) * var(:,9,:) + & q(10,8) * var(:,10,:) + q(11,8) * var(:,11,:) + & q(12,8) * var(:,12,:) ) * idel jl = 9 end if if ( bb(2) == 0 ) then jr = nj - gsize else dvar(:,nj,:) = - ( q(1,1) * var(:,nj,:) + q(2,1) * var(:,nj-1,:) + & q(3,1) * var(:,nj-2,:) + q(4,1) * var(:,nj-3,:) + & q(5,1) * var(:,nj-4,:) + q(6,1) * var(:,nj-5,:) + & q(7,1) * var(:,nj-6,:) + & q(8,1) * var(:,nj-7,:) ) * idel dvar(:,nj-1,:) = - ( q(1,2) * var(:,nj,:) + q(3,2) * var(:,nj-2,:) + & q(4,2) * var(:,nj-3,:) + q(5,2) * var(:,nj-4,:) + & q(6,2) * var(:,nj-5,:) + q(7,2) * var(:,nj-6,:) + & q(8,2) * var(:,nj-7,:) ) * idel dvar(:,nj-2,:) = - ( q(1,3) * var(:,nj,:) + q(2,3) * var(:,nj-1,:) + & q(4,3) * var(:,nj-3,:) + q(5,3) * var(:,nj-4,:) + & q(6,3) * var(:,nj-5,:) + q(7,3) * var(:,nj-6,:) + & q(8,3) * var(:,nj-7,:) ) * idel dvar(:,nj-3,:) = - ( q(1,4) * var(:,nj,:) + q(2,4) * var(:,nj-1,:) + & q(3,4) * var(:,nj-2,:) + q(5,4) * var(:,nj-4,:) + & q(6,4) * var(:,nj-5,:) + q(7,4) * var(:,nj-6,:) + & q(8,4) * var(:,nj-7,:) ) * idel dvar(:,nj-4,:) = - ( q(1,5) * var(:,nj,:) + q(2,5) * var(:,nj-1,:) + & q(3,5) * var(:,nj-2,:) + q(4,5) * var(:,nj-3,:) + & q(6,5) * var(:,nj-5,:) + q(7,5) * var(:,nj-6,:) + & q(8,5) * var(:,nj-7,:) + & q(9,5) * var(:,nj-8,:) ) * idel dvar(:,nj-5,:) = - ( q(1,6) * var(:,nj,:) + q(2,6) * var(:,nj-1,:) + & q(3,6) * var(:,nj-2,:) + q(4,6) * var(:,nj-3,:) + & q(5,6) * var(:,nj-4,:) + q(7,6) * var(:,nj-6,:) + & q(8,6) * var(:,nj-7,:) + q(9,6) * var(:,nj-8,:) + & q(10,6) * var(:,nj-9,:) ) * idel dvar(:,nj-6,:) = - ( q(1,7) * var(:,nj,:) + q(2,7) * var(:,nj-1,:) + & q(3,7) * var(:,nj-2,:) + q(4,7) * var(:,nj-3,:) + & q(5,7) * var(:,nj-4,:) + q(6,7) * var(:,nj-5,:) + & q(8,7) * var(:,nj-7,:) + q(9,7) * var(:,nj-8,:) + & q(10,7) * var(:,nj-9,:) + & q(11,7) * var(:,nj-10,:) ) * idel dvar(:,nj-7,:) = - ( q(1,8) * var(:,nj,:) + q(2,8) * var(:,nj-1,:) + & q(3,8) * var(:,nj-2,:) + q(4,8) * var(:,nj-3,:) + & q(5,8) * var(:,nj-4,:) + q(6,8) * var(:,nj-5,:) + & q(7,8) * var(:,nj-6,:) + q(9,8) * var(:,nj-8,:) + & q(10,8) * var(:,nj-9,:) + & q(11,8) * var(:,nj-10,:) + & q(12,8) * var(:,nj-11,:) ) * idel jr = nj - 8 end if dvar(:,jl:jr,:) = ( a(1) * ( var(:,jl+1:jr+1,:) - & var(:,jl-1:jr-1,:) ) + & a(2) * ( var(:,jl+2:jr+2,:) - & var(:,jl-2:jr-2,:) ) + & a(3) * ( var(:,jl+3:jr+3,:) - & var(:,jl-3:jr-3,:) ) + & a(4) * ( var(:,jl+4:jr+4,:) - & var(:,jl-4:jr-4,:) ) ) * idel end if case (2) direction if ( zero_derivs_z /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then kl = 1 + gsize else dvar(:,:,1) = ( q(1,1) * var(:,:,1) + q(2,1) * var(:,:,2) + & q(3,1) * var(:,:,3) + q(4,1) * var(:,:,4) + & q(5,1) * var(:,:,5) + q(6,1) * var(:,:,6) + & q(7,1) * var(:,:,7) + q(8,1) * var(:,:,8) ) * idel dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(3,2) * var(:,:,3) + & q(4,2) * var(:,:,4) + q(5,2) * var(:,:,5) + & q(6,2) * var(:,:,6) + q(7,2) * var(:,:,7) + & q(8,2) * var(:,:,8) ) * idel dvar(:,:,3) = ( q(1,3) * var(:,:,1) + q(2,3) * var(:,:,2) + & q(4,3) * var(:,:,4) + q(5,3) * var(:,:,5) + & q(6,3) * var(:,:,6) + q(7,3) * var(:,:,7) + & q(8,3) * var(:,:,8) ) * idel dvar(:,:,4) = ( q(1,4) * var(:,:,1) + q(2,4) * var(:,:,2) + & q(3,4) * var(:,:,3) + q(5,4) * var(:,:,5) + & q(6,4) * var(:,:,6) + q(7,4) * var(:,:,7) + & q(8,4) * var(:,:,8) ) * idel dvar(:,:,5) = ( q(1,5) * var(:,:,1) + q(2,5) * var(:,:,2) + & q(3,5) * var(:,:,3) + q(4,5) * var(:,:,4) + & q(6,5) * var(:,:,6) + q(7,5) * var(:,:,7) + & q(8,5) * var(:,:,8) + q(9,5) * var(:,:,9) ) * idel dvar(:,:,6) = ( q(1,6) * var(:,:,1) + q(2,6) * var(:,:,2) + & q(3,6) * var(:,:,3) + q(4,6) * var(:,:,4) + & q(5,6) * var(:,:,5) + q(7,6) * var(:,:,7) + & q(8,6) * var(:,:,8) + q(9,6) * var(:,:,9) + & q(10,6) * var(:,:,10) ) * idel dvar(:,:,7) = ( q(1,7) * var(:,:,1) + q(2,7) * var(:,:,2) + & q(3,7) * var(:,:,3) + q(4,7) * var(:,:,4) + & q(5,7) * var(:,:,5) + q(6,7) * var(:,:,6) + & q(8,7) * var(:,:,8) + q(9,7) * var(:,:,9) + & q(10,7) * var(:,:,10) + q(11,7) * var(:,:,11) ) * idel dvar(:,:,8) = ( q(1,8) * var(:,:,1) + q(2,8) * var(:,:,2) + & q(3,8) * var(:,:,3) + q(4,8) * var(:,:,4) + & q(5,8) * var(:,:,5) + q(6,8) * var(:,:,6) + & q(7,8) * var(:,:,7) + q(9,8) * var(:,:,9) + & q(10,8) * var(:,:,10) + q(11,8) * var(:,:,11) + & q(12,8) * var(:,:,12) ) * idel kl = 9 end if if ( bb(2) == 0 ) then kr = nk - gsize else dvar(:,:,nk) = - ( q(1,1) * var(:,:,nk) + q(2,1) * var(:,:,nk-1) + & q(3,1) * var(:,:,nk-2) + q(4,1) * var(:,:,nk-3) + & q(5,1) * var(:,:,nk-4) + q(6,1) * var(:,:,nk-5) + & q(7,1) * var(:,:,nk-6) + & q(8,1) * var(:,:,nk-7) ) * idel dvar(:,:,nk-1) = - ( q(1,2) * var(:,:,nk) + q(3,2) * var(:,:,nk-2) + & q(4,2) * var(:,:,nk-3) + q(5,2) * var(:,:,nk-4) + & q(6,2) * var(:,:,nk-5) + q(7,2) * var(:,:,nk-6) + & q(8,2) * var(:,:,nk-7) ) * idel dvar(:,:,nk-2) = - ( q(1,3) * var(:,:,nk) + q(2,3) * var(:,:,nk-1) + & q(4,3) * var(:,:,nk-3) + q(5,3) * var(:,:,nk-4) + & q(6,3) * var(:,:,nk-5) + q(7,3) * var(:,:,nk-6) + & q(8,3) * var(:,:,nk-7) ) * idel dvar(:,:,nk-3) = - ( q(1,4) * var(:,:,nk) + q(2,4) * var(:,:,nk-1) + & q(3,4) * var(:,:,nk-2) + q(5,4) * var(:,:,nk-4) + & q(6,4) * var(:,:,nk-5) + q(7,4) * var(:,:,nk-6) + & q(8,4) * var(:,:,nk-7) ) * idel dvar(:,:,nk-4) = - ( q(1,5) * var(:,:,nk) + q(2,5) * var(:,:,nk-1) + & q(3,5) * var(:,:,nk-2) + q(4,5) * var(:,:,nk-3) + & q(6,5) * var(:,:,nk-5) + q(7,5) * var(:,:,nk-6) + & q(8,5) * var(:,:,nk-7) + & q(9,5) * var(:,:,nk-8) ) * idel dvar(:,:,nk-5) = - ( q(1,6) * var(:,:,nk) + q(2,6) * var(:,:,nk-1) + & q(3,6) * var(:,:,nk-2) + q(4,6) * var(:,:,nk-3) + & q(5,6) * var(:,:,nk-4) + q(7,6) * var(:,:,nk-6) + & q(8,6) * var(:,:,nk-7) + q(9,6) * var(:,:,nk-8) + & q(10,6) * var(:,:,nk-9) ) * idel dvar(:,:,nk-6) = - ( q(1,7) * var(:,:,nk) + q(2,7) * var(:,:,nk-1) + & q(3,7) * var(:,:,nk-2) + q(4,7) * var(:,:,nk-3) + & q(5,7) * var(:,:,nk-4) + q(6,7) * var(:,:,nk-5) + & q(8,7) * var(:,:,nk-7) + q(9,7) * var(:,:,nk-8) + & q(10,7) * var(:,:,nk-9) + & q(11,7) * var(:,:,nk-10) ) * idel dvar(:,:,nk-7) = - ( q(1,8) * var(:,:,nk) + q(2,8) * var(:,:,nk-1) + & q(3,8) * var(:,:,nk-2) + q(4,8) * var(:,:,nk-3) + & q(5,8) * var(:,:,nk-4) + q(6,8) * var(:,:,nk-5) + & q(7,8) * var(:,:,nk-6) + q(9,8) * var(:,:,nk-8) + & q(10,8) * var(:,:,nk-9) + & q(11,8) * var(:,:,nk-10) + & q(12,8) * var(:,:,nk-11) ) * idel kr = nk - 8 end if dvar(:,:,kl:kr) = ( a(1) * ( var(:,:,kl+1:kr+1) - & var(:,:,kl-1:kr-1) ) + & a(2) * ( var(:,:,kl+2:kr+2) - & var(:,:,kl-2:kr-2) ) + & a(3) * ( var(:,:,kl+3:kr+3) - & var(:,:,kl-3:kr-3) ) + & a(4) * ( var(:,:,kl+4:kr+4) - & var(:,:,kl-4:kr-4) ) ) * idel end if end select direction end subroutine deriv_gf_8_4_opt subroutine up_deriv_gf_8_4_opt ( var, ni, nj, nk, dir, bb, gsize, delta, up, dvar ) implicit none DECLARE_CCTK_FUNCTIONS DECLARE_CCTK_PARAMETERS CCTK_REAL, parameter :: zero = 0.0 integer, parameter :: wp = kind(zero) CCTK_INT, intent(IN) :: ni, nj, nk CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: var CCTK_INT, intent(IN) :: dir CCTK_INT, intent(IN) :: bb(2) CCTK_INT, intent(IN) :: gsize CCTK_REAL, intent(IN) :: delta CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: up CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar CCTK_REAL, dimension(-4:4), save :: a1, a2 CCTK_REAL, dimension(12,8), save :: q1, q2 CCTK_REAL :: idel CCTK_INT :: il, ir, jl, jr, kl, kr, j, k logical, save :: first = .true. if ( first ) then a1(-4) = 0.007142857142857142857142857142857142857143_wp a1(-3) = -0.06666666666666666666666666666666666666667_wp a1(-2) = 0.3000000000000000000000000000000000000000_wp a1(-1) = -1.000000000000000000000000000000000000000_wp a1(0) = 0.2500000000000000000000000000000000000000_wp a1(1) = 0.6000000000000000000000000000000000000000_wp a1(2) = -0.1000000000000000000000000000000000000000_wp a1(3) = 0.009523809523809523809523809523809523809524_wp a1(4) = 0.0_wp q1(1,1) = -1.683432578685956376544499542432311020539_wp q1(2,1) = 2.206194911574014373726323572790957192479_wp q1(3,1) = 0.01048240012328238642042347285056590509052_wp q1(4,1) = -0.7964177278055455600674641925035752355813_wp q1(5,1) = 0.04375559918421074726330123748705042737925_wp q1(6,1) = 0.3015760468426855559502835132688146245819_wp q1(7,1) = -0.05430139927462465206917597295687512552949_wp q1(8,1) = -0.02785725195806647467919208850462676788091_wp q1(9,1) = 0.0_wp q1(10,1) = 0.0_wp q1(11,1) = 0.0_wp q1(12,1) = 0.0_wp q1(1,2) = -0.4451389901846271204230470438663315357999_wp q1(2,2) = 0.03979384217746298000464189244717138065595_wp q1(3,2) = 0.04117629055092890361160914968898310705568_wp q1(4,2) = 0.5200531283081606472663885367960872642342_wp q1(5,2) = -0.04132093916257534380117029102933732094903_wp q1(6,2) = -0.1571076966553020594965523100929511511562_wp q1(7,2) = 0.02929412303865500620649355002594213960776_wp q1(8,2) = 0.01325024192729698663163651603043611635147_wp q1(9,2) = 0.0_wp q1(10,2) = 0.0_wp q1(11,2) = 0.0_wp q1(12,2) = 0.0_wp q1(1,3) = 0.1544592597396873173299548977341719367844_wp q1(2,3) = -1.020860666961436408391054234781341175626_wp q1(3,3) = 0.7352247001589518809305910119936419247628_wp q1(4,3) = -0.1443377614671085870438942853092495421459_wp q1(5,3) = 0.4029099592852542668089969628980513126845_wp q1(6,3) = -0.1256847607545176255327746274725689708913_wp q1(7,3) = -0.01218362768224537680217842793138395523727_wp q1(8,3) = 0.01047289768141453270035870286867846966855_wp q1(9,3) = 0.0_wp q1(10,3) = 0.0_wp q1(11,3) = 0.0_wp q1(12,3) = 0.0_wp q1(1,4) = 0.1147228529148158318034738218073197976023_wp q1(2,4) = -0.3538780254601914958115698874206531961307_wp q1(3,4) = -0.1858994785495101740456777336107169678288_wp q1(4,4) = 0.1370483808799853486477767415913728669524_wp q1(5,4) = 0.03972858582710281564215102156521588367881_wp q1(6,4) = 0.3165092134256211709442618752193670178495_wp q1(7,4) = -0.05069213274750316093259254157360772183295_wp q1(8,4) = -0.01753939629032033624782329757829768029059_wp q1(9,4) = 0.0_wp q1(10,4) = 0.0_wp q1(11,4) = 0.0_wp q1(12,4) = 0.0_wp q1(1,5) = -0.01395722952056041916292538026746562805620_wp q1(2,5) = 0.01429909551156929498323180030121632055909_wp q1(3,5) = 0.2332633686362482363893930312775209255800_wp q1(4,5) = -1.142300242635165935993053851544552454320_wp q1(5,5) = 0.6057550738330768177827040634066377989296_wp q1(6,5) = 0.2825179352075121763892365873707254541218_wp q1(7,5) = 0.03248906901967492013975174489420235738435_wp q1(8,5) = -0.01206707005235509052833799543828477419916_wp q1(9,5) = 0.0_wp q1(10,5) = 0.0_wp q1(11,5) = 0.0_wp q1(12,5) = 0.0_wp q1(1,6) = -0.06956044896551484587496684695381321009857_wp q1(2,6) = 0.1930764794969705526456194423505107491362_wp q1(3,6) = -0.01938619631973248083178041998353003089641_wp q1(4,6) = -0.2887164590533642022788579933885334942677_wp q1(5,6) = -0.4040701147451219074780049288053335780831_wp q1(6,6) = 0.1955440022503079624043376625896000853565_wp q1(7,6) = 0.4088327186474708718613662440475307694516_wp q1(8,6) = -0.02316927663483720615835459462174938908831_wp q1(9,6) = 0.007449295323821255710641434765318098489772_wp q1(10,6) = 0.0_wp q1(11,6) = 0.0_wp q1(12,6) = 0.0_wp q1(1,7) = 0.01734328280101122760645584447264624051821_wp q1(2,7) = -0.04840773491893482968357096782164797435370_wp q1(3,7) = 0.01113356043999974094579726795398111697087_wp q1(4,7) = 0.03683258325116868830333840353178374785556_wp q1(5,7) = 0.2020929196538139300412937632456998385290_wp q1(6,7) = -0.9993401510385675084229179036326626049945_wp q1(7,7) = 0.2707691940348400441134114620253070105404_wp q1(8,7) = 0.6075690064749917706800934212435751620821_wp q1(9,7) = -0.1083076776139360176453645848101228042161_wp q1(10,7) = 0.01031501691561295406146329379144026706820_wp q1(11,7) = 0.0_wp q1(12,7) = 0.0_wp q1(1,8) = 0.008138876723484902470545040666116789759842_wp q1(2,8) = -0.02002921745130282774392822052620050919550_wp q1(3,8) = -0.002671343688148396036480262140129972226295_wp q1(4,8) = 0.03832298452325717814060194729385919629895_wp q1(5,8) = -0.05168030333817646678101782026786038350804_wp q1(6,8) = 0.2274981280312773565363597229555297756596_wp q1(7,8) = -0.9520792032950265164134838698189786312337_wp q1(8,8) = 0.2476881136110866797237534905735256806167_wp q1(9,8) = 0.5944514726666080313370083773764616334802_wp q1(10,8) = -0.09907524544443467188950139622941027224670_wp q1(11,8) = 0.009435737661374730656142990117086692594923_wp q1(12,8) = 0.0_wp a2(-4) = 0.0_wp a2(-3) = -0.009523809523809523809523809523809523809524_wp a2(-2) = 0.1000000000000000000000000000000000000000_wp a2(-1) = -0.6000000000000000000000000000000000000000_wp a2(0) = -0.2500000000000000000000000000000000000000_wp a2(1) = 1.000000000000000000000000000000000000000_wp a2(2) = -0.3000000000000000000000000000000000000000_wp a2(3) = 0.06666666666666666666666666666666666666667_wp a2(4) = -0.007142857142857142857142857142857142857143_wp q2(1,1) = -1.707654630177840640955211766064430603569_wp q2(2,1) = 2.303083117541551431369172467319435524597_wp q2(3,1) = -0.1348499088280232000438498689421515930856_wp q2(4,1) = -0.6995295218380085024246152979750969034639_wp q2(5,1) = 0.01953354769232648285258901385493084434991_wp q2(6,1) = 0.3015760468426855559502835132688146245819_wp q2(7,1) = -0.05430139927462465206917597295687512552949_wp q2(8,1) = -0.02785725195806647467919208850462676788091_wp q2(9,1) = 0.0_wp q2(10,1) = 0.0_wp q2(11,1) = 0.0_wp q2(12,1) = 0.0_wp q2(1,2) = -0.4264124762187621886561567415382508860795_wp q2(2,2) = -0.03979384217746298000464189244717138065595_wp q2(3,2) = 0.1722618883119834259798412659855476550988_wp q2(4,2) = 0.4170573014959035225484918739916436907718_wp q2(5,2) = -0.003867911230845480267389686373176021508129_wp q2(6,2) = -0.1617893251467682924382748856749713135863_wp q2(7,2) = 0.02929412303865500620649355002594213960776_wp q2(8,2) = 0.01325024192729698663163651603043611635147_wp q2(9,2) = 0.0_wp q2(10,2) = 0.0_wp q2(11,2) = 0.0_wp q2(12,2) = 0.0_wp q2(1,3) = -0.01200671010762254250338646347193566882227_wp q2(2,3) = -0.2440194743406570625021278824861723494614_wp q2(3,3) = -0.7352247001589518809305910119936419247628_wp q2(4,3) = 1.298367310542910198178397511810349706445_wp q2(5,3) = -0.3739312333355250790799293893971175134799_wp q2(6,3) = 0.09626986570856218757834718746890783658430_wp q2(7,3) = -0.03992795599013035344106865479906855617172_wp q2(8,3) = 0.01047289768141453270035870286867846966855_wp q2(9,3) = 0.0_wp q2(10,3) = 0.0_wp q2(11,3) = 0.0_wp q2(12,3) = 0.0_wp q2(1,4) = 0.1306125202632199301974189512671891155098_wp q2(2,4) = -0.4412711958764140369782680994499344446221_wp q2(3,4) = 0.02066619697974310507560894936758416496903_wp q2(4,4) = -0.1370483808799853486477767415913728669524_wp q2(5,4) = 0.2621839287047601931573828340033863343841_wp q2(6,4) = 0.2052815419867924821866459690002817924969_wp q2(7,4) = -0.01891279805069496414470228265386908601791_wp q2(8,4) = -0.02151181312742136084630957994326500976747_wp q2(9,4) = 0.0_wp q2(10,4) = 0.0_wp q2(11,4) = 0.0_wp q2(12,4) = 0.0_wp q2(1,5) = -0.03126451734436261395671692493622670802562_wp q2(2,5) = 0.1527573981019868533335641576513049603144_wp q2(3,5) = -0.2513406904302132178367702194477893135637_wp q2(4,5) = -0.1730921245022430275407273500939319760321_wp q2(5,5) = -0.6057550738330768177827040634066377989296_wp q2(6,5) = 1.251726053340435084841563088821345932409_wp q2(7,5) = -0.4521149900467865340864115058311078817594_wp q2(8,5) = 0.1263912325380624678219943619118038655562_wp q2(9,5) = -0.01730728782380219479379154466876107996942_wp q2(10,5) = 0.0_wp q2(11,5) = 0.0_wp q2(12,5) = 0.0_wp q2(1,6) = -0.06956044896551484587496684695381321009857_wp q2(2,6) = 0.1874895080041046108626383662765221752689_wp q2(3,6) = 0.02530957562319505343206818860837856004222_wp q2(4,6) = -0.4451516608536105722023281234602135625529_wp q2(5,6) = -0.09119971114462916763106466866197344151272_wp q2(6,6) = -0.1955440022503079624043376625896000853565_wp q2(7,6) = 0.7217031222479636117083065041908909060220_wp q2(8,6) = -0.1796044784350835760818247246934294573735_wp q2(9,6) = 0.05214506726674878997449004335722668942840_wp q2(10,6) = -0.005586971492865941782981076073988573867329_wp q2(11,6) = 0.0_wp q2(12,6) = 0.0_wp q2(1,7) = 0.01734328280101122760645584447264624051821_wp q2(2,7) = -0.04840773491893482968357096782164797435370_wp q2(3,7) = 0.003397297753290025399699797610400916669720_wp q2(4,7) = 0.09872268474484641267211816628042535026479_wp q2(5,7) = -0.01452243557405810524943540637454576990329_wp q2(6,7) = -0.5661094405828234378414595643921713881299_wp q2(7,7) = -0.2707691940348400441134114620253070105404_wp q2(8,7) = 1.040799716930735841261551760484066378947_wp q2(9,7) = -0.3249230328418080529360937544303684126484_wp q2(10,7) = 0.07220511840929067843024305654008186947743_wp q2(11,7) = -0.007736262686709715546097470343580200301153_wp q2(12,7) = 0.0_wp q2(1,8) = 0.008138876723484902470545040666116789759842_wp q2(2,8) = -0.02002921745130282774392822052620050919550_wp q2(3,8) = -0.002671343688148396036480262140129972226295_wp q2(4,8) = 0.03124618127722613014849470470604417685276_wp q2(5,8) = 0.004934122630071917155840120434659772061498_wp q2(6,8) = 0.02934763714240801275735693049670923116616_wp q2(7,8) = -0.5557782215172878288554782849013375422469_wp q2(8,8) = -0.2476881136110866797237534905735256806167_wp q2(9,8) = 0.9907524544443467188950139622941027224670_wp q2(10,8) = -0.2972257363333040156685041886882308167401_wp q2(11,8) = 0.06605016362962311459300093081960684816446_wp q2(12,8) = -0.007076803246031047992107242587815019446193_wp first = .false. end if idel = 1.0_wp / delta direction: select case (dir) case (0) direction if ( bb(1) == 0 ) then il = 1 + gsize else where ( up(1,:,:) < zero ) dvar(1,:,:) = ( q1(1,1) * var(1,:,:) + q1(2,1) * var(2,:,:) + & q1(3,1) * var(3,:,:) + q1(4,1) * var(4,:,:) + & q1(5,1) * var(5,:,:) + q1(6,1) * var(6,:,:) + & q1(7,1) * var(7,:,:) + q1(8,1) * var(8,:,:) ) * idel elsewhere dvar(1,:,:) = ( q2(1,1) * var(1,:,:) + q2(2,1) * var(2,:,:) + & q2(3,1) * var(3,:,:) + q2(4,1) * var(4,:,:) + & q2(5,1) * var(5,:,:) + q2(6,1) * var(6,:,:) + & q2(7,1) * var(7,:,:) + q2(8,1) * var(8,:,:) ) * idel end where where ( up(2,:,:) < zero ) dvar(2,:,:) = ( q1(1,2) * var(1,:,:) + q1(2,2) * var(2,:,:) + & q1(3,2) * var(3,:,:) + q1(4,2) * var(4,:,:) + & q1(5,2) * var(5,:,:) + q1(6,2) * var(6,:,:) + & q1(7,2) * var(7,:,:) + q1(8,2) * var(8,:,:) ) * idel elsewhere dvar(2,:,:) = ( q2(1,2) * var(1,:,:) + q2(2,2) * var(2,:,:) + & q2(3,2) * var(3,:,:) + q2(4,2) * var(4,:,:) + & q2(5,2) * var(5,:,:) + q2(6,2) * var(6,:,:) + & q2(7,2) * var(7,:,:) + q2(8,2) * var(8,:,:) ) * idel end where where ( up(3,:,:) < zero ) dvar(3,:,:) = ( q1(1,3) * var(1,:,:) + q1(2,3) * var(2,:,:) + & q1(3,3) * var(3,:,:) + q1(4,3) * var(4,:,:) + & q1(5,3) * var(5,:,:) + q1(6,3) * var(6,:,:) + & q1(7,3) * var(7,:,:) + q1(8,3) * var(8,:,:) ) * idel elsewhere dvar(3,:,:) = ( q2(1,3) * var(1,:,:) + q2(2,3) * var(2,:,:) + & q2(3,3) * var(3,:,:) + q2(4,3) * var(4,:,:) + & q2(5,3) * var(5,:,:) + q2(6,3) * var(6,:,:) + & q2(7,3) * var(7,:,:) + q2(8,3) * var(8,:,:) ) * idel end where where ( up(4,:,:) < zero ) dvar(4,:,:) = ( q1(1,4) * var(1,:,:) + q1(2,4) * var(2,:,:) + & q1(3,4) * var(3,:,:) + q1(4,4) * var(4,:,:) + & q1(5,4) * var(5,:,:) + q1(6,4) * var(6,:,:) + & q1(7,4) * var(7,:,:) + q1(8,4) * var(8,:,:) ) * idel elsewhere dvar(4,:,:) = ( q2(1,4) * var(1,:,:) + q2(2,4) * var(2,:,:) + & q2(3,4) * var(3,:,:) + q2(4,4) * var(4,:,:) + & q2(5,4) * var(5,:,:) + q2(6,4) * var(6,:,:) + & q2(7,4) * var(7,:,:) + q2(8,4) * var(8,:,:) ) * idel end where where ( up(5,:,:) < zero ) dvar(5,:,:) = ( q1(1,5) * var(1,:,:) + q1(2,5) * var(2,:,:) + & q1(3,5) * var(3,:,:) + q1(4,5) * var(4,:,:) + & q1(5,5) * var(5,:,:) + q1(6,5) * var(6,:,:) + & q1(7,5) * var(7,:,:) + q1(8,5) * var(8,:,:) + & q1(9,5) * var(9,:,:) ) * idel elsewhere dvar(5,:,:) = ( q2(1,5) * var(1,:,:) + q2(2,5) * var(2,:,:) + & q2(3,5) * var(3,:,:) + q2(4,5) * var(4,:,:) + & q2(5,5) * var(5,:,:) + q2(6,5) * var(6,:,:) + & q2(7,5) * var(7,:,:) + q2(8,5) * var(8,:,:) + & q2(9,5) * var(9,:,:) ) * idel end where where ( up(6,:,:) < zero ) dvar(6,:,:) = ( q1(1,6) * var(1,:,:) + q1(2,6) * var(2,:,:) + & q1(3,6) * var(3,:,:) + q1(4,6) * var(4,:,:) + & q1(5,6) * var(5,:,:) + q1(6,6) * var(6,:,:) + & q1(7,6) * var(7,:,:) + q1(8,6) * var(8,:,:) + & q1(9,6) * var(9,:,:) + q1(10,6) * var(10,:,:) ) * idel elsewhere dvar(6,:,:) = ( q2(1,6) * var(1,:,:) + q2(2,6) * var(2,:,:) + & q2(3,6) * var(3,:,:) + q2(4,6) * var(4,:,:) + & q2(5,6) * var(5,:,:) + q2(6,6) * var(6,:,:) + & q2(7,6) * var(7,:,:) + q2(8,6) * var(8,:,:) + & q2(9,6) * var(9,:,:) + q2(10,6) * var(10,:,:) ) * idel end where where ( up(7,:,:) < zero ) dvar(7,:,:) = ( q1(1,7) * var(1,:,:) + q1(2,7) * var(2,:,:) + & q1(3,7) * var(3,:,:) + q1(4,7) * var(4,:,:) + & q1(5,7) * var(5,:,:) + q1(6,7) * var(6,:,:) + & q1(7,7) * var(7,:,:) + q1(8,7) * var(8,:,:) + & q1(9,7) * var(9,:,:) + q1(10,7) * var(10,:,:) + & q1(11,7) * var(11,:,:) ) * idel elsewhere dvar(7,:,:) = ( q2(1,7) * var(1,:,:) + q2(2,7) * var(2,:,:) + & q2(3,7) * var(3,:,:) + q2(4,7) * var(4,:,:) + & q2(5,7) * var(5,:,:) + q2(6,7) * var(6,:,:) + & q2(7,7) * var(7,:,:) + q2(8,7) * var(8,:,:) + & q2(9,7) * var(9,:,:) + q2(10,7) * var(10,:,:) + & q2(11,7) * var(11,:,:) ) * idel end where where ( up(8,:,:) < zero ) dvar(8,:,:) = ( q1(1,8) * var(1,:,:) + q1(2,8) * var(2,:,:) + & q1(3,8) * var(3,:,:) + q1(4,8) * var(4,:,:) + & q1(5,8) * var(5,:,:) + q1(6,8) * var(6,:,:) + & q1(7,8) * var(7,:,:) + q1(8,8) * var(8,:,:) + & q1(9,8) * var(9,:,:) + q1(10,8) * var(10,:,:) + & q1(11,8) * var(11,:,:) + & q1(12,8) * var(12,:,:) ) * idel elsewhere dvar(8,:,:) = ( q2(1,8) * var(1,:,:) + q2(2,8) * var(2,:,:) + & q2(3,8) * var(3,:,:) + q2(4,8) * var(4,:,:) + & q2(5,8) * var(5,:,:) + q2(6,8) * var(6,:,:) + & q2(7,8) * var(7,:,:) + q2(8,8) * var(8,:,:) + & q2(9,8) * var(9,:,:) + q2(10,8) * var(10,:,:) + & q2(11,8) * var(11,:,:) + & q2(12,8) * var(12,:,:) ) * idel end where il = 9 end if if ( bb(2) == 0 ) then ir = ni - gsize else where ( up(ni,:,:) < zero ) dvar(ni,:,:) = - ( q2(1,1) * var(ni,:,:) + & q2(2,1) * var(ni-1,:,:) + & q2(3,1) * var(ni-2,:,:) + & q2(4,1) * var(ni-3,:,:) + & q2(5,1) * var(ni-4,:,:) + & q2(6,1) * var(ni-5,:,:) + & q2(7,1) * var(ni-6,:,:) + & q2(8,1) * var(ni-7,:,:) ) * idel elsewhere dvar(ni,:,:) = - ( q1(1,1) * var(ni,:,:) + & q1(2,1) * var(ni-1,:,:) + & q1(3,1) * var(ni-2,:,:) + & q1(4,1) * var(ni-3,:,:) + & q1(5,1) * var(ni-4,:,:) + & q1(6,1) * var(ni-5,:,:) + & q1(7,1) * var(ni-6,:,:) + & q1(8,1) * var(ni-7,:,:) ) * idel end where where ( up(ni-1,:,:) < zero ) dvar(ni-1,:,:) = - ( q2(1,2) * var(ni,:,:) + & q2(2,2) * var(ni-1,:,:) + & q2(3,2) * var(ni-2,:,:) + & q2(4,2) * var(ni-3,:,:) + & q2(5,2) * var(ni-4,:,:) + & q2(6,2) * var(ni-5,:,:) + & q2(7,2) * var(ni-6,:,:) + & q2(8,2) * var(ni-7,:,:) ) * idel elsewhere dvar(ni-1,:,:) = - ( q1(1,2) * var(ni,:,:) + & q1(2,2) * var(ni-1,:,:) + & q1(3,2) * var(ni-2,:,:) + & q1(4,2) * var(ni-3,:,:) + & q1(5,2) * var(ni-4,:,:) + & q1(6,2) * var(ni-5,:,:) + & q1(7,2) * var(ni-6,:,:) + & q1(8,2) * var(ni-7,:,:) ) * idel end where where ( up(ni-2,:,:) < zero ) dvar(ni-2,:,:) = - ( q2(1,3) * var(ni,:,:) + & q2(2,3) * var(ni-1,:,:) + & q2(3,3) * var(ni-2,:,:) + & q2(4,3) * var(ni-3,:,:) + & q2(5,3) * var(ni-4,:,:) + & q2(6,3) * var(ni-5,:,:) + & q2(7,3) * var(ni-6,:,:) + & q2(8,3) * var(ni-7,:,:) ) * idel elsewhere dvar(ni-2,:,:) = - ( q1(1,3) * var(ni,:,:) + & q1(2,3) * var(ni-1,:,:) + & q1(3,3) * var(ni-2,:,:) + & q1(4,3) * var(ni-3,:,:) + & q1(5,3) * var(ni-4,:,:) + & q1(6,3) * var(ni-5,:,:) + & q1(7,3) * var(ni-6,:,:) + & q1(8,3) * var(ni-7,:,:) ) * idel end where where ( up(ni-3,:,:) < zero ) dvar(ni-3,:,:) = - ( q2(1,4) * var(ni,:,:) + & q2(2,4) * var(ni-1,:,:) + & q2(3,4) * var(ni-2,:,:) + & q2(4,4) * var(ni-3,:,:) + & q2(5,4) * var(ni-4,:,:) + & q2(6,4) * var(ni-5,:,:) + & q2(7,4) * var(ni-6,:,:) + & q2(8,4) * var(ni-7,:,:) ) * idel elsewhere dvar(ni-3,:,:) = - ( q1(1,4) * var(ni,:,:) + & q1(2,4) * var(ni-1,:,:) + & q1(3,4) * var(ni-2,:,:) + & q1(4,4) * var(ni-3,:,:) + & q1(5,4) * var(ni-4,:,:) + & q1(6,4) * var(ni-5,:,:) + & q1(7,4) * var(ni-6,:,:) + & q1(8,4) * var(ni-7,:,:) ) * idel end where where ( up(ni-4,:,:) < zero ) dvar(ni-4,:,:) = - ( q2(1,5) * var(ni,:,:) + & q2(2,5) * var(ni-1,:,:) + & q2(3,5) * var(ni-2,:,:) + & q2(4,5) * var(ni-3,:,:) + & q2(5,5) * var(ni-4,:,:) + & q2(6,5) * var(ni-5,:,:) + & q2(7,5) * var(ni-6,:,:) + & q2(8,5) * var(ni-7,:,:) + & q2(9,5) * var(ni-8,:,:) ) * idel elsewhere dvar(ni-4,:,:) = - ( q1(1,5) * var(ni,:,:) + & q1(2,5) * var(ni-1,:,:) + & q1(3,5) * var(ni-2,:,:) + & q1(4,5) * var(ni-3,:,:) + & q1(5,5) * var(ni-4,:,:) + & q1(6,5) * var(ni-5,:,:) + & q1(7,5) * var(ni-6,:,:) + & q1(8,5) * var(ni-7,:,:) + & q1(9,5) * var(ni-8,:,:) ) * idel end where where ( up(ni-5,:,:) < zero ) dvar(ni-5,:,:) = - ( q2(1,6) * var(ni,:,:) + & q2(2,6) * var(ni-1,:,:) + & q2(3,6) * var(ni-2,:,:) + & q2(4,6) * var(ni-3,:,:) + & q2(5,6) * var(ni-4,:,:) + & q2(6,6) * var(ni-5,:,:) + & q2(7,6) * var(ni-6,:,:) + & q2(8,6) * var(ni-7,:,:) + & q2(9,6) * var(ni-8,:,:) + & q2(10,6) * var(ni-9,:,:) ) * idel elsewhere dvar(ni-5,:,:) = - ( q1(1,6) * var(ni,:,:) + & q1(2,6) * var(ni-1,:,:) + & q1(3,6) * var(ni-2,:,:) + & q1(4,6) * var(ni-3,:,:) + & q1(5,6) * var(ni-4,:,:) + & q1(6,6) * var(ni-5,:,:) + & q1(7,6) * var(ni-6,:,:) + & q1(8,6) * var(ni-7,:,:) + & q1(9,6) * var(ni-8,:,:) + & q1(10,6) * var(ni-9,:,:) ) * idel end where where ( up(ni-6,:,:) < zero ) dvar(ni-6,:,:) = - ( q2(1,7) * var(ni,:,:) + & q2(2,7) * var(ni-1,:,:) + & q2(3,7) * var(ni-2,:,:) + & q2(4,7) * var(ni-3,:,:) + & q2(5,7) * var(ni-4,:,:) + & q2(6,7) * var(ni-5,:,:) + & q2(7,7) * var(ni-6,:,:) + & q2(8,7) * var(ni-7,:,:) + & q2(9,7) * var(ni-8,:,:) + & q2(10,7) * var(ni-9,:,:) + & q2(11,7) * var(ni-10,:,:) ) * idel elsewhere dvar(ni-6,:,:) = - ( q1(1,7) * var(ni,:,:) + & q1(2,7) * var(ni-1,:,:) + & q1(3,7) * var(ni-2,:,:) + & q1(4,7) * var(ni-3,:,:) + & q1(5,7) * var(ni-4,:,:) + & q1(6,7) * var(ni-5,:,:) + & q1(7,7) * var(ni-6,:,:) + & q1(8,7) * var(ni-7,:,:) + & q1(9,7) * var(ni-8,:,:) + & q1(10,7) * var(ni-9,:,:) + & q1(11,7) * var(ni-10,:,:) ) * idel end where where ( up(ni-7,:,:) < zero ) dvar(ni-7,:,:) = - ( q2(1,8) * var(ni,:,:) + & q2(2,8) * var(ni-1,:,:) + & q2(3,8) * var(ni-2,:,:) + & q2(4,8) * var(ni-3,:,:) + & q2(5,8) * var(ni-4,:,:) + & q2(6,8) * var(ni-5,:,:) + & q2(7,8) * var(ni-6,:,:) + & q2(8,8) * var(ni-7,:,:) + & q2(9,8) * var(ni-8,:,:) + & q2(10,8) * var(ni-9,:,:) + & q2(11,8) * var(ni-10,:,:) + & q2(12,8) * var(ni-11,:,:) ) * idel elsewhere dvar(ni-7,:,:) = - ( q1(1,8) * var(ni,:,:) + & q1(2,8) * var(ni-1,:,:) + & q1(3,8) * var(ni-2,:,:) + & q1(4,8) * var(ni-3,:,:) + & q1(5,8) * var(ni-4,:,:) + & q1(6,8) * var(ni-5,:,:) + & q1(7,8) * var(ni-6,:,:) + & q1(8,8) * var(ni-7,:,:) + & q1(9,8) * var(ni-8,:,:) + & q1(10,8) * var(ni-9,:,:) + & q1(11,8) * var(ni-10,:,:) + & q1(12,8) * var(ni-11,:,:) ) * idel end where ir = ni - 8 end if where ( up(il:ir,:,:) < zero ) dvar(il:ir,:,:) = ( a1(-4) * var(il-4:ir-4,:,:) + & a1(-3) * var(il-3:ir-3,:,:) + & a1(-2) * var(il-2:ir-2,:,:) + & a1(-1) * var(il-1:ir-1,:,:) + & a1(0) * var(il:ir,:,:) + & a1(1) * var(il+1:ir+1,:,:) + & a1(2) * var(il+2:ir+2,:,:) + & a1(3) * var(il+3:ir+3,:,:) + & a1(4) * var(il+4:ir+4,:,:) ) * idel elsewhere dvar(il:ir,:,:) = ( a2(-4) * var(il-4:ir-4,:,:) + & a2(-3) * var(il-3:ir-3,:,:) + & a2(-2) * var(il-2:ir-2,:,:) + & a2(-1) * var(il-1:ir-1,:,:) + & a2(0) * var(il:ir,:,:) + & a2(1) * var(il+1:ir+1,:,:) + & a2(2) * var(il+2:ir+2,:,:) + & a2(3) * var(il+3:ir+3,:,:) + & a2(4) * var(il+4:ir+4,:,:) ) * idel end where case (1) direction if ( zero_derivs_y /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then jl = 1 + gsize else where ( up(:,1,:) < zero ) dvar(:,1,:) = ( q1(1,1) * var(:,1,:) + q1(2,1) * var(:,2,:) + & q1(3,1) * var(:,3,:) + q1(4,1) * var(:,4,:) + & q1(5,1) * var(:,5,:) + q1(6,1) * var(:,6,:) + & q1(7,1) * var(:,7,:) + q1(8,1) * var(:,8,:) ) * idel elsewhere dvar(:,1,:) = ( q2(1,1) * var(:,1,:) + q2(2,1) * var(:,2,:) + & q2(3,1) * var(:,3,:) + q2(4,1) * var(:,4,:) + & q2(5,1) * var(:,5,:) + q2(6,1) * var(:,6,:) + & q2(7,1) * var(:,7,:) + q2(8,1) * var(:,8,:) ) * idel end where where ( up(:,2,:) < zero ) dvar(:,2,:) = ( q1(1,2) * var(:,1,:) + q1(2,2) * var(:,2,:) + & q1(3,2) * var(:,3,:) + q1(4,2) * var(:,4,:) + & q1(5,2) * var(:,5,:) + q1(6,2) * var(:,6,:) + & q1(7,2) * var(:,7,:) + q1(8,2) * var(:,8,:) ) * idel elsewhere dvar(:,2,:) = ( q2(1,2) * var(:,1,:) + q2(2,2) * var(:,2,:) + & q2(3,2) * var(:,3,:) + q2(4,2) * var(:,4,:) + & q2(5,2) * var(:,5,:) + q2(6,2) * var(:,6,:) + & q2(7,2) * var(:,7,:) + q2(8,2) * var(:,8,:) ) * idel end where where ( up(:,3,:) < zero ) dvar(:,3,:) = ( q1(1,3) * var(:,1,:) + q1(2,3) * var(:,2,:) + & q1(3,3) * var(:,3,:) + q1(4,3) * var(:,4,:) + & q1(5,3) * var(:,5,:) + q1(6,3) * var(:,6,:) + & q1(7,3) * var(:,7,:) + q1(8,3) * var(:,8,:) ) * idel elsewhere dvar(:,3,:) = ( q2(1,3) * var(:,1,:) + q2(2,3) * var(:,2,:) + & q2(3,3) * var(:,3,:) + q2(4,3) * var(:,4,:) + & q2(5,3) * var(:,5,:) + q2(6,3) * var(:,6,:) + & q2(7,3) * var(:,7,:) + q2(8,3) * var(:,8,:) ) * idel end where where ( up(:,4,:) < zero ) dvar(:,4,:) = ( q1(1,4) * var(:,1,:) + q1(2,4) * var(:,2,:) + & q1(3,4) * var(:,3,:) + q1(4,4) * var(:,4,:) + & q1(5,4) * var(:,5,:) + q1(6,4) * var(:,6,:) + & q1(7,4) * var(:,7,:) + q1(8,4) * var(:,8,:) ) * idel elsewhere dvar(:,4,:) = ( q2(1,4) * var(:,1,:) + q2(2,4) * var(:,2,:) + & q2(3,4) * var(:,3,:) + q2(4,4) * var(:,4,:) + & q2(5,4) * var(:,5,:) + q2(6,4) * var(:,6,:) + & q2(7,4) * var(:,7,:) + q2(8,4) * var(:,8,:) ) * idel end where where ( up(:,5,:) < zero ) dvar(:,5,:) = ( q1(1,5) * var(:,1,:) + q1(2,5) * var(:,2,:) + & q1(3,5) * var(:,3,:) + q1(4,5) * var(:,4,:) + & q1(5,5) * var(:,5,:) + q1(6,5) * var(:,6,:) + & q1(7,5) * var(:,7,:) + q1(8,5) * var(:,8,:) + & q1(9,5) * var(:,9,:) ) * idel elsewhere dvar(:,5,:) = ( q2(1,5) * var(:,1,:) + q2(2,5) * var(:,2,:) + & q2(3,5) * var(:,3,:) + q2(4,5) * var(:,4,:) + & q2(5,5) * var(:,5,:) + q2(6,5) * var(:,6,:) + & q2(7,5) * var(:,7,:) + q2(8,5) * var(:,8,:) + & q2(9,5) * var(:,9,:) ) * idel end where where ( up(:,6,:) < zero ) dvar(:,6,:) = ( q1(1,6) * var(:,1,:) + q1(2,6) * var(:,2,:) + & q1(3,6) * var(:,3,:) + q1(4,6) * var(:,4,:) + & q1(5,6) * var(:,5,:) + q1(6,6) * var(:,6,:) + & q1(7,6) * var(:,7,:) + q1(8,6) * var(:,8,:) + & q1(9,6) * var(:,9,:) + q1(10,6) * var(:,10,:) ) * idel elsewhere dvar(:,6,:) = ( q2(1,6) * var(:,1,:) + q2(2,6) * var(:,2,:) + & q2(3,6) * var(:,3,:) + q2(4,6) * var(:,4,:) + & q2(5,6) * var(:,5,:) + q2(6,6) * var(:,6,:) + & q2(7,6) * var(:,7,:) + q2(8,6) * var(:,8,:) + & q2(9,6) * var(:,9,:) + q2(10,6) * var(:,10,:) ) * idel end where where ( up(:,7,:) < zero ) dvar(:,7,:) = ( q1(1,7) * var(:,1,:) + q1(2,7) * var(:,2,:) + & q1(3,7) * var(:,3,:) + q1(4,7) * var(:,4,:) + & q1(5,7) * var(:,5,:) + q1(6,7) * var(:,6,:) + & q1(7,7) * var(:,7,:) + q1(8,7) * var(:,8,:) + & q1(9,7) * var(:,9,:) + q1(10,7) * var(:,10,:) + & q1(11,7) * var(:,11,:) ) * idel elsewhere dvar(:,7,:) = ( q2(1,7) * var(:,1,:) + q2(2,7) * var(:,2,:) + & q2(3,7) * var(:,3,:) + q2(4,7) * var(:,4,:) + & q2(5,7) * var(:,5,:) + q2(6,7) * var(:,6,:) + & q2(7,7) * var(:,7,:) + q2(8,7) * var(:,8,:) + & q2(9,7) * var(:,9,:) + q2(10,7) * var(:,10,:) + & q2(11,7) * var(:,11,:) ) * idel end where where ( up(:,8,:) < zero ) dvar(:,8,:) = ( q1(1,8) * var(:,1,:) + q1(2,8) * var(:,2,:) + & q1(3,8) * var(:,3,:) + q1(4,8) * var(:,4,:) + & q1(5,8) * var(:,5,:) + q1(6,8) * var(:,6,:) + & q1(7,8) * var(:,7,:) + q1(8,8) * var(:,8,:) + & q1(9,8) * var(:,9,:) + q1(10,8) * var(:,10,:) + & q1(11,8) * var(:,11,:) + & q1(12,8) * var(:,12,:) ) * idel elsewhere dvar(:,8,:) = ( q2(1,8) * var(:,1,:) + q2(2,8) * var(:,2,:) + & q2(3,8) * var(:,3,:) + q2(4,8) * var(:,4,:) + & q2(5,8) * var(:,5,:) + q2(6,8) * var(:,6,:) + & q2(7,8) * var(:,7,:) + q2(8,8) * var(:,8,:) + & q2(9,8) * var(:,9,:) + q2(10,8) * var(:,10,:) + & q2(11,8) * var(:,11,:) + & q2(12,8) * var(:,12,:) ) * idel end where jl = 9 end if if ( bb(2) == 0 ) then jr = nj - gsize else where ( up(:,nj,:) < zero ) dvar(:,nj,:) = - ( q2(1,1) * var(:,nj,:) + & q2(2,1) * var(:,nj-1,:) + & q2(3,1) * var(:,nj-2,:) + & q2(4,1) * var(:,nj-3,:) + & q2(5,1) * var(:,nj-4,:) + & q2(6,1) * var(:,nj-5,:) + & q2(7,1) * var(:,nj-6,:) + & q2(8,1) * var(:,nj-7,:) ) * idel elsewhere dvar(:,nj,:) = - ( q1(1,1) * var(:,nj,:) + & q1(2,1) * var(:,nj-1,:) + & q1(3,1) * var(:,nj-2,:) + & q1(4,1) * var(:,nj-3,:) + & q1(5,1) * var(:,nj-4,:) + & q1(6,1) * var(:,nj-5,:) + & q1(7,1) * var(:,nj-6,:) + & q1(8,1) * var(:,nj-7,:) ) * idel end where where ( up(:,nj-1,:) < zero ) dvar(:,nj-1,:) = - ( q2(1,2) * var(:,nj,:) + & q2(2,2) * var(:,nj-1,:) + & q2(3,2) * var(:,nj-2,:) + & q2(4,2) * var(:,nj-3,:) + & q2(5,2) * var(:,nj-4,:) + & q2(6,2) * var(:,nj-5,:) + & q2(7,2) * var(:,nj-6,:) + & q2(8,2) * var(:,nj-7,:) ) * idel elsewhere dvar(:,nj-1,:) = - ( q1(1,2) * var(:,nj,:) + & q1(2,2) * var(:,nj-1,:) + & q1(3,2) * var(:,nj-2,:) + & q1(4,2) * var(:,nj-3,:) + & q1(5,2) * var(:,nj-4,:) + & q1(6,2) * var(:,nj-5,:) + & q1(7,2) * var(:,nj-6,:) + & q1(8,2) * var(:,nj-7,:) ) * idel end where where ( up(:,nj-2,:) < zero ) dvar(:,nj-2,:) = - ( q2(1,3) * var(:,nj,:) + & q2(2,3) * var(:,nj-1,:) + & q2(3,3) * var(:,nj-2,:) + & q2(4,3) * var(:,nj-3,:) + & q2(5,3) * var(:,nj-4,:) + & q2(6,3) * var(:,nj-5,:) + & q2(7,3) * var(:,nj-6,:) + & q2(8,3) * var(:,nj-7,:) ) * idel elsewhere dvar(:,nj-2,:) = - ( q1(1,3) * var(:,nj,:) + & q1(2,3) * var(:,nj-1,:) + & q1(3,3) * var(:,nj-2,:) + & q1(4,3) * var(:,nj-3,:) + & q1(5,3) * var(:,nj-4,:) + & q1(6,3) * var(:,nj-5,:) + & q1(7,3) * var(:,nj-6,:) + & q1(8,3) * var(:,nj-7,:) ) * idel end where where ( up(:,nj-3,:) < zero ) dvar(:,nj-3,:) = - ( q2(1,4) * var(:,nj,:) + & q2(2,4) * var(:,nj-1,:) + & q2(3,4) * var(:,nj-2,:) + & q2(4,4) * var(:,nj-3,:) + & q2(5,4) * var(:,nj-4,:) + & q2(6,4) * var(:,nj-5,:) + & q2(7,4) * var(:,nj-6,:) + & q2(8,4) * var(:,nj-7,:) ) * idel elsewhere dvar(:,nj-3,:) = - ( q1(1,4) * var(:,nj,:) + & q1(2,4) * var(:,nj-1,:) + & q1(3,4) * var(:,nj-2,:) + & q1(4,4) * var(:,nj-3,:) + & q1(5,4) * var(:,nj-4,:) + & q1(6,4) * var(:,nj-5,:) + & q1(7,4) * var(:,nj-6,:) + & q1(8,4) * var(:,nj-7,:) ) * idel end where where ( up(:,nj-4,:) < zero ) dvar(:,nj-4,:) = - ( q2(1,5) * var(:,nj,:) + & q2(2,5) * var(:,nj-1,:) + & q2(3,5) * var(:,nj-2,:) + & q2(4,5) * var(:,nj-3,:) + & q2(5,5) * var(:,nj-4,:) + & q2(6,5) * var(:,nj-5,:) + & q2(7,5) * var(:,nj-6,:) + & q2(8,5) * var(:,nj-7,:) + & q2(9,5) * var(:,nj-8,:) ) * idel elsewhere dvar(:,nj-4,:) = - ( q1(1,5) * var(:,nj,:) + & q1(2,5) * var(:,nj-1,:) + & q1(3,5) * var(:,nj-2,:) + & q1(4,5) * var(:,nj-3,:) + & q1(5,5) * var(:,nj-4,:) + & q1(6,5) * var(:,nj-5,:) + & q1(7,5) * var(:,nj-6,:) + & q1(8,5) * var(:,nj-7,:) + & q1(9,5) * var(:,nj-8,:) ) * idel end where where ( up(:,nj-5,:) < zero ) dvar(:,nj-5,:) = - ( q2(1,6) * var(:,nj,:) + & q2(2,6) * var(:,nj-1,:) + & q2(3,6) * var(:,nj-2,:) + & q2(4,6) * var(:,nj-3,:) + & q2(5,6) * var(:,nj-4,:) + & q2(6,6) * var(:,nj-5,:) + & q2(7,6) * var(:,nj-6,:) + & q2(8,6) * var(:,nj-7,:) + & q2(9,6) * var(:,nj-8,:) + & q2(10,6) * var(:,nj-9,:) ) * idel elsewhere dvar(:,nj-5,:) = - ( q1(1,6) * var(:,nj,:) + & q1(2,6) * var(:,nj-1,:) + & q1(3,6) * var(:,nj-2,:) + & q1(4,6) * var(:,nj-3,:) + & q1(5,6) * var(:,nj-4,:) + & q1(6,6) * var(:,nj-5,:) + & q1(7,6) * var(:,nj-6,:) + & q1(8,6) * var(:,nj-7,:) + & q1(9,6) * var(:,nj-8,:) + & q1(10,6) * var(:,nj-9,:) ) * idel end where where ( up(:,nj-6,:) < zero ) dvar(:,nj-6,:) = - ( q2(1,7) * var(:,nj,:) + & q2(2,7) * var(:,nj-1,:) + & q2(3,7) * var(:,nj-2,:) + & q2(4,7) * var(:,nj-3,:) + & q2(5,7) * var(:,nj-4,:) + & q2(6,7) * var(:,nj-5,:) + & q2(7,7) * var(:,nj-6,:) + & q2(8,7) * var(:,nj-7,:) + & q2(9,7) * var(:,nj-8,:) + & q2(10,7) * var(:,nj-9,:) + & q2(11,7) * var(:,nj-10,:) ) * idel elsewhere dvar(:,nj-6,:) = - ( q1(1,7) * var(:,nj,:) + & q1(2,7) * var(:,nj-1,:) + & q1(3,7) * var(:,nj-2,:) + & q1(4,7) * var(:,nj-3,:) + & q1(5,7) * var(:,nj-4,:) + & q1(6,7) * var(:,nj-5,:) + & q1(7,7) * var(:,nj-6,:) + & q1(8,7) * var(:,nj-7,:) + & q1(9,7) * var(:,nj-8,:) + & q1(10,7) * var(:,nj-9,:) + & q1(11,7) * var(:,nj-10,:) ) * idel end where where ( up(:,nj-7,:) < zero ) dvar(:,nj-7,:) = - ( q2(1,8) * var(:,nj,:) + & q2(2,8) * var(:,nj-1,:) + & q2(3,8) * var(:,nj-2,:) + & q2(4,8) * var(:,nj-3,:) + & q2(5,8) * var(:,nj-4,:) + & q2(6,8) * var(:,nj-5,:) + & q2(7,8) * var(:,nj-6,:) + & q2(8,8) * var(:,nj-7,:) + & q2(9,8) * var(:,nj-8,:) + & q2(10,8) * var(:,nj-9,:) + & q2(11,8) * var(:,nj-10,:) + & q2(12,8) * var(:,nj-11,:) ) * idel elsewhere dvar(:,nj-7,:) = - ( q1(1,8) * var(:,nj,:) + & q1(2,8) * var(:,nj-1,:) + & q1(3,8) * var(:,nj-2,:) + & q1(4,8) * var(:,nj-3,:) + & q1(5,8) * var(:,nj-4,:) + & q1(6,8) * var(:,nj-5,:) + & q1(7,8) * var(:,nj-6,:) + & q1(8,8) * var(:,nj-7,:) + & q1(9,8) * var(:,nj-8,:) + & q1(10,8) * var(:,nj-9,:) + & q1(11,8) * var(:,nj-10,:) + & q1(12,8) * var(:,nj-11,:) ) * idel end where jr = nj - 8 end if where ( up(:,jl:jr,:) < zero ) dvar(:,jl:jr,:) = ( a1(-4) * var(:,jl-4:jr-4,:) + & a1(-3) * var(:,jl-3:jr-3,:) + & a1(-2) * var(:,jl-2:jr-2,:) + & a1(-1) * var(:,jl-1:jr-1,:) + & a1(0) * var(:,jl:jr,:) + & a1(1) * var(:,jl+1:jr+1,:) + & a1(2) * var(:,jl+2:jr+2,:) + & a1(3) * var(:,jl+3:jr+3,:) + & a1(4) * var(:,jl+4:jr+4,:) ) * idel elsewhere dvar(:,jl:jr,:) = ( a2(-4) * var(:,jl-4:jr-4,:) + & a2(-3) * var(:,jl-3:jr-3,:) + & a2(-2) * var(:,jl-2:jr-2,:) + & a2(-1) * var(:,jl-1:jr-1,:) + & a2(0) * var(:,jl:jr,:) + & a2(1) * var(:,jl+1:jr+1,:) + & a2(2) * var(:,jl+2:jr+2,:) + & a2(3) * var(:,jl+3:jr+3,:) + & a2(4) * var(:,jl+4:jr+4,:) ) * idel end where end if case (2) direction if ( zero_derivs_z /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then kl = 1 + gsize else where ( up(:,:,1) < zero ) dvar(:,:,1) = ( q2(1,1) * var(:,:,1) + q2(2,1) * var(:,:,2) + & q2(3,1) * var(:,:,3) + q2(4,1) * var(:,:,4) + & q2(5,1) * var(:,:,5) + q2(6,1) * var(:,:,6) + & q2(7,1) * var(:,:,7) + q2(8,1) * var(:,:,8) ) * idel elsewhere dvar(:,:,1) = ( q2(1,1) * var(:,:,1) + q2(2,1) * var(:,:,2) + & q2(3,1) * var(:,:,3) + q2(4,1) * var(:,:,4) + & q2(5,1) * var(:,:,5) + q2(6,1) * var(:,:,6) + & q2(7,1) * var(:,:,7) + q2(8,1) * var(:,:,8) ) * idel end where where ( up(:,:,2) < zero ) dvar(:,:,2) = ( q1(1,2) * var(:,:,1) + q1(2,2) * var(:,:,2) + & q1(3,2) * var(:,:,3) + q1(4,2) * var(:,:,4) + & q1(5,2) * var(:,:,5) + q1(6,2) * var(:,:,6) + & q1(7,2) * var(:,:,7) + q1(8,2) * var(:,:,8) ) * idel elsewhere dvar(:,:,2) = ( q2(1,2) * var(:,:,1) + q2(2,2) * var(:,:,2) + & q2(3,2) * var(:,:,3) + q2(4,2) * var(:,:,4) + & q2(5,2) * var(:,:,5) + q2(6,2) * var(:,:,6) + & q2(7,2) * var(:,:,7) + q2(8,2) * var(:,:,8) ) * idel end where where ( up(:,:,3) < zero ) dvar(:,:,3) = ( q1(1,3) * var(:,:,1) + q1(2,3) * var(:,:,2) + & q1(3,3) * var(:,:,3) + q1(4,3) * var(:,:,4) + & q1(5,3) * var(:,:,5) + q1(6,3) * var(:,:,6) + & q1(7,3) * var(:,:,7) + q1(8,3) * var(:,:,8) ) * idel elsewhere dvar(:,:,3) = ( q2(1,3) * var(:,:,1) + q2(2,3) * var(:,:,2) + & q2(3,3) * var(:,:,3) + q2(4,3) * var(:,:,4) + & q2(5,3) * var(:,:,5) + q2(6,3) * var(:,:,6) + & q2(7,3) * var(:,:,7) + q2(8,3) * var(:,:,8) ) * idel end where where ( up(:,:,4) < zero ) dvar(:,:,4) = ( q1(1,4) * var(:,:,1) + q1(2,4) * var(:,:,2) + & q1(3,4) * var(:,:,3) + q1(4,4) * var(:,:,4) + & q1(5,4) * var(:,:,5) + q1(6,4) * var(:,:,6) + & q1(7,4) * var(:,:,7) + q1(8,4) * var(:,:,8) ) * idel elsewhere dvar(:,:,4) = ( q2(1,4) * var(:,:,1) + q2(2,4) * var(:,:,2) + & q2(3,4) * var(:,:,3) + q2(4,4) * var(:,:,4) + & q2(5,4) * var(:,:,5) + q2(6,4) * var(:,:,6) + & q2(7,4) * var(:,:,7) + q2(8,4) * var(:,:,8) ) * idel end where where ( up(:,:,5) < zero ) dvar(:,:,5) = ( q1(1,5) * var(:,:,1) + q1(2,5) * var(:,:,2) + & q1(3,5) * var(:,:,3) + q1(4,5) * var(:,:,4) + & q1(5,5) * var(:,:,5) + q1(6,5) * var(:,:,6) + & q1(7,5) * var(:,:,7) + q1(8,5) * var(:,:,8) + & q1(9,5) * var(:,:,9) ) * idel elsewhere dvar(:,:,5) = ( q2(1,5) * var(:,:,1) + q2(2,5) * var(:,:,2) + & q2(3,5) * var(:,:,3) + q2(4,5) * var(:,:,4) + & q2(5,5) * var(:,:,5) + q2(6,5) * var(:,:,6) + & q2(7,5) * var(:,:,7) + q2(8,5) * var(:,:,8) + & q2(9,5) * var(:,:,9) ) * idel end where where ( up(:,:,6) < zero ) dvar(:,:,6) = ( q1(1,6) * var(:,:,1) + q1(2,6) * var(:,:,2) + & q1(3,6) * var(:,:,3) + q1(4,6) * var(:,:,4) + & q1(5,6) * var(:,:,5) + q1(6,6) * var(:,:,6) + & q1(7,6) * var(:,:,7) + q1(8,6) * var(:,:,8) + & q1(9,6) * var(:,:,9) + q1(10,6) * var(:,:,10) ) * idel elsewhere dvar(:,:,6) = ( q2(1,6) * var(:,:,1) + q2(2,6) * var(:,:,2) + & q2(3,6) * var(:,:,3) + q2(4,6) * var(:,:,4) + & q2(5,6) * var(:,:,5) + q2(6,6) * var(:,:,6) + & q2(7,6) * var(:,:,7) + q2(8,6) * var(:,:,8) + & q2(9,6) * var(:,:,9) + q2(10,6) * var(:,:,10) ) * idel end where where ( up(:,:,7) < zero ) dvar(:,:,7) = ( q1(1,7) * var(:,:,1) + q1(2,7) * var(:,:,2) + & q1(3,7) * var(:,:,3) + q1(4,7) * var(:,:,4) + & q1(5,7) * var(:,:,5) + q1(6,7) * var(:,:,6) + & q1(7,7) * var(:,:,7) + q1(8,7) * var(:,:,8) + & q1(9,7) * var(:,:,9) + q1(10,7) * var(:,:,10) + & q1(11,7) * var(:,:,11) ) * idel elsewhere dvar(:,:,7) = ( q2(1,7) * var(:,:,1) + q2(2,7) * var(:,:,2) + & q2(3,7) * var(:,:,3) + q2(4,7) * var(:,:,4) + & q2(5,7) * var(:,:,5) + q2(6,7) * var(:,:,6) + & q2(7,7) * var(:,:,7) + q2(8,7) * var(:,:,8) + & q2(9,7) * var(:,:,9) + q2(10,7) * var(:,:,10) + & q2(11,7) * var(:,:,11) ) * idel end where where ( up(:,:,8) < zero ) dvar(:,:,8) = ( q1(1,8) * var(:,:,1) + q1(2,8) * var(:,:,2) + & q1(3,8) * var(:,:,3) + q1(4,8) * var(:,:,4) + & q1(5,8) * var(:,:,5) + q1(6,8) * var(:,:,6) + & q1(7,8) * var(:,:,7) + q1(8,8) * var(:,:,8) + & q1(9,8) * var(:,:,9) + q1(10,8) * var(:,:,10) + & q1(11,8) * var(:,:,11) + & q1(12,8) * var(:,:,12) ) * idel elsewhere dvar(:,:,8) = ( q2(1,8) * var(:,:,1) + q2(2,8) * var(:,:,2) + & q2(3,8) * var(:,:,3) + q2(4,8) * var(:,:,4) + & q2(5,8) * var(:,:,5) + q2(6,8) * var(:,:,6) + & q2(7,8) * var(:,:,7) + q2(8,8) * var(:,:,8) + & q2(9,8) * var(:,:,9) + q2(10,8) * var(:,:,10) + & q2(11,8) * var(:,:,11) + & q2(12,8) * var(:,:,12) ) * idel end where kl = 9 end if if ( bb(2) == 0 ) then kr = nk - gsize else where ( up(:,:,nk) < zero ) dvar(:,:,nk) = - ( q2(1,1) * var(:,:,nk) + & q2(2,1) * var(:,:,nk-1) + & q2(3,1) * var(:,:,nk-2) + & q2(4,1) * var(:,:,nk-3) + & q2(5,1) * var(:,:,nk-4) + & q2(6,1) * var(:,:,nk-5) + & q2(7,1) * var(:,:,nk-6) + & q2(8,1) * var(:,:,nk-7) ) * idel elsewhere dvar(:,:,nk) = - ( q1(1,1) * var(:,:,nk) + & q1(2,1) * var(:,:,nk-1) + & q1(3,1) * var(:,:,nk-2) + & q1(4,1) * var(:,:,nk-3) + & q1(5,1) * var(:,:,nk-4) + & q1(6,1) * var(:,:,nk-5) + & q1(7,1) * var(:,:,nk-6) + & q1(8,1) * var(:,:,nk-7) ) * idel end where where ( up(:,:,nk-1) < zero ) dvar(:,:,nk-1) = - ( q2(1,2) * var(:,:,nk) + & q2(2,2) * var(:,:,nk-1) + & q2(3,2) * var(:,:,nk-2) + & q2(4,2) * var(:,:,nk-3) + & q2(5,2) * var(:,:,nk-4) + & q2(6,2) * var(:,:,nk-5) + & q2(7,2) * var(:,:,nk-6) + & q2(8,2) * var(:,:,nk-7) ) * idel elsewhere dvar(:,:,nk-1) = - ( q1(1,2) * var(:,:,nk) + & q1(2,2) * var(:,:,nk-1) + & q1(3,2) * var(:,:,nk-2) + & q1(4,2) * var(:,:,nk-3) + & q1(5,2) * var(:,:,nk-4) + & q1(6,2) * var(:,:,nk-5) + & q1(7,2) * var(:,:,nk-6) + & q1(8,2) * var(:,:,nk-7) ) * idel end where where ( up(:,:,nk-2) < zero ) dvar(:,:,nk-2) = - ( q2(1,3) * var(:,:,nk) + & q2(2,3) * var(:,:,nk-1) + & q2(3,3) * var(:,:,nk-2) + & q2(4,3) * var(:,:,nk-3) + & q2(5,3) * var(:,:,nk-4) + & q2(6,3) * var(:,:,nk-5) + & q2(7,3) * var(:,:,nk-6) + & q2(8,3) * var(:,:,nk-7) ) * idel elsewhere dvar(:,:,nk-2) = - ( q1(1,3) * var(:,:,nk) + & q1(2,3) * var(:,:,nk-1) + & q1(3,3) * var(:,:,nk-2) + & q1(4,3) * var(:,:,nk-3) + & q1(5,3) * var(:,:,nk-4) + & q1(6,3) * var(:,:,nk-5) + & q1(7,3) * var(:,:,nk-6) + & q1(8,3) * var(:,:,nk-7) ) * idel end where where ( up(:,:,nk-3) < zero ) dvar(:,:,nk-3) = - ( q2(1,4) * var(:,:,nk) + & q2(2,4) * var(:,:,nk-1) + & q2(3,4) * var(:,:,nk-2) + & q2(4,4) * var(:,:,nk-3) + & q2(5,4) * var(:,:,nk-4) + & q2(6,4) * var(:,:,nk-5) + & q2(7,4) * var(:,:,nk-6) + & q2(8,4) * var(:,:,nk-7) ) * idel elsewhere dvar(:,:,nk-3) = - ( q1(1,4) * var(:,:,nk) + & q1(2,4) * var(:,:,nk-1) + & q1(3,4) * var(:,:,nk-2) + & q1(4,4) * var(:,:,nk-3) + & q1(5,4) * var(:,:,nk-4) + & q1(6,4) * var(:,:,nk-5) + & q1(7,4) * var(:,:,nk-6) + & q1(8,4) * var(:,:,nk-7) ) * idel end where where ( up(:,:,nk-4) < zero ) dvar(:,:,nk-4) = - ( q2(1,5) * var(:,:,nk) + & q2(2,5) * var(:,:,nk-1) + & q2(3,5) * var(:,:,nk-2) + & q2(4,5) * var(:,:,nk-3) + & q2(5,5) * var(:,:,nk-4) + & q2(6,5) * var(:,:,nk-5) + & q2(7,5) * var(:,:,nk-6) + & q2(8,5) * var(:,:,nk-7) + & q2(9,5) * var(:,:,nk-8) ) * idel elsewhere dvar(:,:,nk-4) = - ( q1(1,5) * var(:,:,nk) + & q1(2,5) * var(:,:,nk-1) + & q1(3,5) * var(:,:,nk-2) + & q1(4,5) * var(:,:,nk-3) + & q1(5,5) * var(:,:,nk-4) + & q1(6,5) * var(:,:,nk-5) + & q1(7,5) * var(:,:,nk-6) + & q1(8,5) * var(:,:,nk-7) + & q1(9,5) * var(:,:,nk-8) ) * idel end where where ( up(:,:,nk-5) < zero ) dvar(:,:,nk-5) = - ( q2(1,6) * var(:,:,nk) + & q2(2,6) * var(:,:,nk-1) + & q2(3,6) * var(:,:,nk-2) + & q2(4,6) * var(:,:,nk-3) + & q2(5,6) * var(:,:,nk-4) + & q2(6,6) * var(:,:,nk-5) + & q2(7,6) * var(:,:,nk-6) + & q2(8,6) * var(:,:,nk-7) + & q2(9,6) * var(:,:,nk-8) + & q2(10,6) * var(:,:,nk-9) ) * idel elsewhere dvar(:,:,nk-5) = - ( q1(1,6) * var(:,:,nk) + & q1(2,6) * var(:,:,nk-1) + & q1(3,6) * var(:,:,nk-2) + & q1(4,6) * var(:,:,nk-3) + & q1(5,6) * var(:,:,nk-4) + & q1(6,6) * var(:,:,nk-5) + & q1(7,6) * var(:,:,nk-6) + & q1(8,6) * var(:,:,nk-7) + & q1(9,6) * var(:,:,nk-8) + & q1(10,6) * var(:,:,nk-9) ) * idel end where where ( up(:,:,nk-6) < zero ) dvar(:,:,nk-6) = - ( q2(1,7) * var(:,:,nk) + & q2(2,7) * var(:,:,nk-1) + & q2(3,7) * var(:,:,nk-2) + & q2(4,7) * var(:,:,nk-3) + & q2(5,7) * var(:,:,nk-4) + & q2(6,7) * var(:,:,nk-5) + & q2(7,7) * var(:,:,nk-6) + & q2(8,7) * var(:,:,nk-7) + & q2(9,7) * var(:,:,nk-8) + & q2(10,7) * var(:,:,nk-9) + & q2(11,7) * var(:,:,nk-10) ) * idel elsewhere dvar(:,:,nk-6) = - ( q1(1,7) * var(:,:,nk) + & q1(2,7) * var(:,:,nk-1) + & q1(3,7) * var(:,:,nk-2) + & q1(4,7) * var(:,:,nk-3) + & q1(5,7) * var(:,:,nk-4) + & q1(6,7) * var(:,:,nk-5) + & q1(7,7) * var(:,:,nk-6) + & q1(8,7) * var(:,:,nk-7) + & q1(9,7) * var(:,:,nk-8) + & q1(10,7) * var(:,:,nk-9) + & q1(11,7) * var(:,:,nk-10) ) * idel end where where ( up(:,:,nk-7) < zero ) dvar(:,:,nk-7) = - ( q2(1,8) * var(:,:,nk) + & q2(2,8) * var(:,:,nk-1) + & q2(3,8) * var(:,:,nk-2) + & q2(4,8) * var(:,:,nk-3) + & q2(5,8) * var(:,:,nk-4) + & q2(6,8) * var(:,:,nk-5) + & q2(7,8) * var(:,:,nk-6) + & q2(8,8) * var(:,:,nk-7) + & q2(9,8) * var(:,:,nk-8) + & q2(10,8) * var(:,:,nk-9) + & q2(11,8) * var(:,:,nk-10) + & q2(12,8) * var(:,:,nk-11) ) * idel elsewhere dvar(:,:,nk-7) = - ( q1(1,8) * var(:,:,nk) + & q1(2,8) * var(:,:,nk-1) + & q1(3,8) * var(:,:,nk-2) + & q1(4,8) * var(:,:,nk-3) + & q1(5,8) * var(:,:,nk-4) + & q1(6,8) * var(:,:,nk-5) + & q1(7,8) * var(:,:,nk-6) + & q1(8,8) * var(:,:,nk-7) + & q1(9,8) * var(:,:,nk-8) + & q1(10,8) * var(:,:,nk-9) + & q1(11,8) * var(:,:,nk-10) + & q1(12,8) * var(:,:,nk-11) ) * idel end where kr = nk - 8 end if where ( up(:,:,kl:kr) < zero ) dvar(:,:,kl:kr) = ( a1(-4) * var(:,:,kl-4:kr-4) + & a1(-3) * var(:,:,kl-3:kr-3) + & a1(-2) * var(:,:,kl-2:kr-2) + & a1(-1) * var(:,:,kl-1:kr-1) + & a1(0) * var(:,:,kl:kr) + & a1(1) * var(:,:,kl+1:kr+1) + & a1(2) * var(:,:,kl+2:kr+2) + & a1(3) * var(:,:,kl+3:kr+3) + & a1(4) * var(:,:,kl+4:kr+4) ) * idel elsewhere dvar(:,:,kl:kr) = ( a2(-4) * var(:,:,kl-4:kr-4) + & a2(-3) * var(:,:,kl-3:kr-3) + & a2(-2) * var(:,:,kl-2:kr-2) + & a2(-1) * var(:,:,kl-1:kr-1) + & a2(0) * var(:,:,kl:kr) + & a2(1) * var(:,:,kl+1:kr+1) + & a2(2) * var(:,:,kl+2:kr+2) + & a2(3) * var(:,:,kl+3:kr+3) + & a2(4) * var(:,:,kl+4:kr+4) ) * idel end where end if end select direction end subroutine up_deriv_gf_8_4_opt