1*59599516SKenneth E. Jansen!!***************************************************** 2*59599516SKenneth E. Jansen! 3*59599516SKenneth E. Jansen! Turbo AMG (ramg) interface functions 4*59599516SKenneth E. Jansen! 1. V_cycle 5*59599516SKenneth E. Jansen! 2. plain 6*59599516SKenneth E. Jansen! 3. driver 7*59599516SKenneth E. Jansen! 8*59599516SKenneth E. Jansen!***************************************************** 9*59599516SKenneth E. Jansen 10*59599516SKenneth E. Jansen !******************************************* 11*59599516SKenneth E. Jansen ! ramg_V_cycle: 12*59599516SKenneth E. Jansen ! One V_cycle call. 13*59599516SKenneth E. Jansen ! Now use SAMG, later will use own code 14*59599516SKenneth E. Jansen !******************************************* 15*59599516SKenneth E. Jansen recursive subroutine ramg_V_cycle( 16*59599516SKenneth E. Jansen & ramg_sol,ramg_rhs,clevel, 17*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper 18*59599516SKenneth E. Jansen & )!ramg_V_cycle 19*59599516SKenneth E. Jansen use ramg_data 20*59599516SKenneth E. Jansen include "common.h" 21*59599516SKenneth E. Jansen include "mpif.h" 22*59599516SKenneth E. Jansen 23*59599516SKenneth E. Jansen !Variable Declaration 24*59599516SKenneth E. Jansen ! arguments 25*59599516SKenneth E. Jansen integer, intent(in) :: clevel 26*59599516SKenneth E. Jansen real(kind=8),intent(inout),dimension(amg_nshg(clevel)) 27*59599516SKenneth E. Jansen & :: ramg_sol,ramg_rhs 28*59599516SKenneth E. Jansen 29*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg+1) :: colm 30*59599516SKenneth E. Jansen integer,intent(in),dimension(nnz_tot) :: rowp 31*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK 32*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP 33*59599516SKenneth E. Jansen integer,intent(in),dimension(nlwork) :: ilwork 34*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg) :: iBC,iper 35*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC 36*59599516SKenneth E. Jansen 37*59599516SKenneth E. Jansen 38*59599516SKenneth E. Jansen ! local 39*59599516SKenneth E. Jansen real(kind=8),dimension(amg_nshg(clevel)) :: myvF,myvE 40*59599516SKenneth E. Jansen real(kind=8),dimension(:),allocatable :: myvC,myvCS 41*59599516SKenneth E. Jansen real(kind=8) :: cpusec(10) 42*59599516SKenneth E. Jansen 43*59599516SKenneth E. Jansen ! In scale 44*59599516SKenneth E. Jansen 45*59599516SKenneth E. Jansen if (clevel.ne.1) then 46*59599516SKenneth E. Jansen ramg_rhs = ramg_rhs*amg_ppeDiag(clevel)%p 47*59599516SKenneth E. Jansen endif 48*59599516SKenneth E. Jansen 49*59599516SKenneth E. Jansen if (ramg_levelx-clevel .eq. 0) then 50*59599516SKenneth E. Jansen ! finest level solver 51*59599516SKenneth E. Jansen call ramg_direct_solve(ramg_rhs,ramg_sol, 52*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP, 53*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper) 54*59599516SKenneth E. Jansen !********************** 55*59599516SKenneth E. Jansen ! higher level smoother 56*59599516SKenneth E. Jansen !********************** 57*59599516SKenneth E. Jansen else if (ramg_levelx-clevel .gt. 0) then 58*59599516SKenneth E. Jansen allocate(myvC(amg_nshg(clevel+1))) 59*59599516SKenneth E. Jansen allocate(myvCS(amg_nshg(clevel+1))) 60*59599516SKenneth E. Jansen 61*59599516SKenneth E. Jansen ! smoothing (pre) 62*59599516SKenneth E. Jansen myvF = 0 ! initial guess is zero 63*59599516SKenneth E. Jansen call ramg_smoother(clevel,ramg_sol,myvF,ramg_rhs, 64*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP, 65*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper,2) !dx1 66*59599516SKenneth E. Jansen 67*59599516SKenneth E. Jansen ! restriction 68*59599516SKenneth E. Jansen call ramg_calcAv_g(clevel,myvF,ramg_sol,colm,rowp,lhsK,lhsP, 69*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper,1) 70*59599516SKenneth E. Jansen myvF = myvF - ramg_rhs! r2 71*59599516SKenneth E. Jansen call ramg_calcIvFC(clevel,clevel+1,myvF,myvC) 72*59599516SKenneth E. Jansen ! up one level, solve 73*59599516SKenneth E. Jansen myvCS = 0 74*59599516SKenneth E. Jansen call ramg_V_cycle(myvCS,myvC,clevel+1, 75*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 76*59599516SKenneth E. Jansen ! prolongation 77*59599516SKenneth E. Jansen call ramg_calcIvCF(clevel,clevel+1,myvCS,myvF, 78*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper) ! v2hat 79*59599516SKenneth E. Jansen ramg_sol = ramg_sol - myvF 80*59599516SKenneth E. Jansen 81*59599516SKenneth E. Jansen ! smoothing (post) 82*59599516SKenneth E. Jansen call ramg_smoother(clevel,myvF,ramg_sol,ramg_rhs, 83*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP, 84*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper,1) 85*59599516SKenneth E. Jansen ramg_sol = myvF 86*59599516SKenneth E. Jansen 87*59599516SKenneth E. Jansen deallocate(myVC) 88*59599516SKenneth E. Jansen deallocate(myVCS) 89*59599516SKenneth E. Jansen end if 90*59599516SKenneth E. Jansen 91*59599516SKenneth E. Jansen ! Out scale 92*59599516SKenneth E. Jansen if (clevel.ne.1) then 93*59599516SKenneth E. Jansen ramg_sol = ramg_sol * amg_ppeDiag(clevel)%p 94*59599516SKenneth E. Jansen endif 95*59599516SKenneth E. Jansen 96*59599516SKenneth E. Jansen return 97*59599516SKenneth E. Jansen 98*59599516SKenneth E. Jansen end subroutine ramg_V_cycle 99*59599516SKenneth E. Jansen 100*59599516SKenneth E. Jansen 101*59599516SKenneth E. Jansen!***************************************************** 102*59599516SKenneth E. Jansen! ramg Cg solver 103*59599516SKenneth E. Jansen! using ramg_V_cycle as preconditioner 104*59599516SKenneth E. Jansen!***************************************************** 105*59599516SKenneth E. Jansen subroutine ramg_CG( 106*59599516SKenneth E. Jansen & sol,rhs,level, 107*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper,iterNum) 108*59599516SKenneth E. Jansen use ramg_data 109*59599516SKenneth E. Jansen include "common.h" 110*59599516SKenneth E. Jansen include "mpif.h" 111*59599516SKenneth E. Jansen include "auxmpi.h" 112*59599516SKenneth E. Jansen 113*59599516SKenneth E. Jansen integer,intent(in) :: level 114*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(amg_nshg(level)) :: rhs 115*59599516SKenneth E. Jansen real(kind=8),intent(inout),dimension(amg_nshg(level)) :: sol 116*59599516SKenneth E. Jansen 117*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg+1) :: colm 118*59599516SKenneth E. Jansen integer,intent(in),dimension(nnz_tot) :: rowp 119*59599516SKenneth E. Jansen 120*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK 121*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP 122*59599516SKenneth E. Jansen 123*59599516SKenneth E. Jansen integer,intent(in),dimension(nlwork) :: ilwork 124*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg) :: iBC,iper 125*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC 126*59599516SKenneth E. Jansen integer,intent(inout) :: iterNum 127*59599516SKenneth E. Jansen 128*59599516SKenneth E. Jansen ! local 129*59599516SKenneth E. Jansen real(kind=8),dimension(amg_nshg(level)) :: cgP,cgQ,cgZ,cgR 130*59599516SKenneth E. Jansen real(kind=8) :: rz,rz_0,pq,alpha,beta 131*59599516SKenneth E. Jansen real(kind=8) :: tmp,norm_0,norm_1,norm_e,norm_c 132*59599516SKenneth E. Jansen real(kind=8) :: mres_0,mres_n 133*59599516SKenneth E. Jansen 134*59599516SKenneth E. Jansen cgR = rhs 135*59599516SKenneth E. Jansen call ramg_L2_norm_p(level,cgR,1,norm_0) 136*59599516SKenneth E. Jansen norm_c = norm_0 137*59599516SKenneth E. Jansen norm_e = norm_0*1e-7 138*59599516SKenneth E. Jansen 139*59599516SKenneth E. Jansen cgZ = cgR 140*59599516SKenneth E. Jansen 141*59599516SKenneth E. Jansen call ramg_dot_p(level,cgZ,cgR,1,rz) 142*59599516SKenneth E. Jansen rz_0 = rz 143*59599516SKenneth E. Jansen cgP = cgZ 144*59599516SKenneth E. Jansen sol = 0 145*59599516SKenneth E. Jansen 146*59599516SKenneth E. Jansen iterNum = 1 147*59599516SKenneth E. Jansen 148*59599516SKenneth E. Jansen do while ( (norm_c .gt. norm_e).and.(iterNum.lt.500) ) 149*59599516SKenneth E. Jansen call ramg_calcAv_g(level,cgQ,cgP,colm,rowp,lhsK,lhsP, 150*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper,1) 151*59599516SKenneth E. Jansen call ramg_dot_p(level,cgP,cgQ,1,pq) 152*59599516SKenneth E. Jansen alpha = rz/pq 153*59599516SKenneth E. Jansen sol = sol + alpha*cgP 154*59599516SKenneth E. Jansen cgR = cgR - alpha*cgQ 155*59599516SKenneth E. Jansen call ramg_L2_norm_p(level,cgR,1,tmp) 156*59599516SKenneth E. Jansen norm_c = tmp 157*59599516SKenneth E. Jansen cgZ = cgR 158*59599516SKenneth E. Jansen tmp = rz 159*59599516SKenneth E. Jansen call ramg_dot_p(level,cgZ,cgR,1,rz) 160*59599516SKenneth E. Jansen beta = rz/tmp 161*59599516SKenneth E. Jansen cgP = beta*cgP+cgZ 162*59599516SKenneth E. Jansen iterNum = iterNum + 1 163*59599516SKenneth E. Jansen enddo 164*59599516SKenneth E. Jansen 165*59599516SKenneth E. Jansen end subroutine ! ramg_CG 166*59599516SKenneth E. Jansen 167*59599516SKenneth E. Jansen subroutine ramg_PCG( 168*59599516SKenneth E. Jansen & sol,rhs,run_mode, 169*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper 170*59599516SKenneth E. Jansen & )!ramg_CG 171*59599516SKenneth E. Jansen use ramg_data 172*59599516SKenneth E. Jansen include "common.h" 173*59599516SKenneth E. Jansen include "mpif.h" 174*59599516SKenneth E. Jansen 175*59599516SKenneth E. Jansen real(kind=8),intent(inout),dimension(nshg) :: sol 176*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg) :: rhs 177*59599516SKenneth E. Jansen integer,intent(in) :: run_mode 178*59599516SKenneth E. Jansen 179*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg+1) :: colm 180*59599516SKenneth E. Jansen integer,intent(in),dimension(nnz_tot) :: rowp 181*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK 182*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP 183*59599516SKenneth E. Jansen integer,intent(in),dimension(nlwork) :: ilwork 184*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg) :: iBC,iper 185*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC 186*59599516SKenneth E. Jansen 187*59599516SKenneth E. Jansen 188*59599516SKenneth E. Jansen ! local 189*59599516SKenneth E. Jansen real(kind=8),dimension(amg_nshg(1)) :: cgP,cgQ,cgZ,cgR 190*59599516SKenneth E. Jansen real(kind=8) :: rz,rz_0,pq,alpha,beta 191*59599516SKenneth E. Jansen real(kind=8) :: tmp,norm_0,norm_p,norm_e,norm_c 192*59599516SKenneth E. Jansen real(kind=8) :: mres_0,mres_n 193*59599516SKenneth E. Jansen real(kind=8),dimension(nshg,4) :: diag 194*59599516SKenneth E. Jansen integer :: iterNum 195*59599516SKenneth E. Jansen 196*59599516SKenneth E. Jansen ! controls 197*59599516SKenneth E. Jansen 198*59599516SKenneth E. Jansen real(kind=8),dimension(maxiters) :: rconv 199*59599516SKenneth E. Jansen real(kind=8) :: avgconv 200*59599516SKenneth E. Jansen logical :: vcycle,gcycle,restart 201*59599516SKenneth E. Jansen character cflagv,cflagg 202*59599516SKenneth E. Jansen 203*59599516SKenneth E. Jansen cgR = rhs 204*59599516SKenneth E. Jansen call ramg_L2_norm_p(1,cgR,1,norm_0) 205*59599516SKenneth E. Jansen norm_c = norm_0 206*59599516SKenneth E. Jansen norm_e = norm_0 * epstol(2) 207*59599516SKenneth E. Jansen write(*,*)'norm_0 ', norm_0 208*59599516SKenneth E. Jansen 209*59599516SKenneth E. Jansen vcycle = .true. 210*59599516SKenneth E. Jansen gcycle = .true. 211*59599516SKenneth E. Jansen restart = .false. 212*59599516SKenneth E. Jansen 213*59599516SKenneth E. Jansen ! amg preconditioning control 214*59599516SKenneth E. Jansen 215*59599516SKenneth E. Jansen if (vcycle) then 216*59599516SKenneth E. Jansen call ramg_V_cycle(cgZ,cgR,1, 217*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 218*59599516SKenneth E. Jansen if (gcycle) then 219*59599516SKenneth E. Jansen call ramg_G_cycle(cgZ,cgR,1, 220*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 221*59599516SKenneth E. Jansen endif 222*59599516SKenneth E. Jansen else 223*59599516SKenneth E. Jansen cgZ = cgR 224*59599516SKenneth E. Jansen endif 225*59599516SKenneth E. Jansen 226*59599516SKenneth E. Jansen call ramg_dot_p(1,cgZ,cgR,1,rz) 227*59599516SKenneth E. Jansen rz_0 = rz 228*59599516SKenneth E. Jansen 229*59599516SKenneth E. Jansen cgP = cgZ 230*59599516SKenneth E. Jansen sol = 0 231*59599516SKenneth E. Jansen 232*59599516SKenneth E. Jansen iterNum = 1 233*59599516SKenneth E. Jansen rconv(1) = 1.0 234*59599516SKenneth E. Jansen norm_p=norm_c 235*59599516SKenneth E. Jansen avgconv = 0 236*59599516SKenneth E. Jansen 237*59599516SKenneth E. Jansen do while ( (norm_c.gt.norm_e).and.(iterNum.lt.maxiters)) 238*59599516SKenneth E. Jansen !do while ((iterNum.lt.maxiters)) 239*59599516SKenneth E. Jansen call ramg_calcAv_g(1,cgQ,cgP,colm,rowp,lhsK,lhsP, 240*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper,1) 241*59599516SKenneth E. Jansen call ramg_dot_p(1,cgP,cgQ,1,pq) 242*59599516SKenneth E. Jansen alpha = rz/pq 243*59599516SKenneth E. Jansen sol = sol + alpha*cgP 244*59599516SKenneth E. Jansen cgR = cgR - alpha*cgQ 245*59599516SKenneth E. Jansen call ramg_L2_norm_p(1,cgR,1,tmp) 246*59599516SKenneth E. Jansen norm_c = tmp 247*59599516SKenneth E. Jansen if (iterNum.eq.1) rconv(1) = norm_p/norm_c 248*59599516SKenneth E. Jansen 249*59599516SKenneth E. Jansen ! amg preconditioning control 250*59599516SKenneth E. Jansen if (vcycle) then 251*59599516SKenneth E. Jansen call ramg_V_cycle(cgZ,cgR,1, 252*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 253*59599516SKenneth E. Jansen if (gcycle) then 254*59599516SKenneth E. Jansen call ramg_G_cycle(cgZ,cgR,1, 255*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 256*59599516SKenneth E. Jansen endif 257*59599516SKenneth E. Jansen else 258*59599516SKenneth E. Jansen cgZ = cgR 259*59599516SKenneth E. Jansen endif 260*59599516SKenneth E. Jansen 261*59599516SKenneth E. Jansen tmp = rz 262*59599516SKenneth E. Jansen call ramg_dot_p(1,cgZ,cgR,1,rz) 263*59599516SKenneth E. Jansen beta = rz/tmp 264*59599516SKenneth E. Jansen cgP = beta*cgP + cgZ 265*59599516SKenneth E. Jansen iterNum = iterNum + 1 266*59599516SKenneth E. Jansen rconv(iterNum) = norm_p/norm_c 267*59599516SKenneth E. Jansen norm_p = norm_c 268*59599516SKenneth E. Jansen cflagv='n' 269*59599516SKenneth E. Jansen cflagg='n' 270*59599516SKenneth E. Jansen if (vcycle) cflagv='v' 271*59599516SKenneth E. Jansen if (gcycle) cflagg='g' 272*59599516SKenneth E. Jansen if (myrank.eq.master) then 273*59599516SKenneth E. Jansen write(*,"((A7)(I4)(E14.4)(T30)(F8.4)(T40)(A1)(A1))") 274*59599516SKenneth E. Jansen & 'AMGCG: ', 275*59599516SKenneth E. Jansen & iterNum,norm_c/norm_0,rconv(iterNum),cflagv,cflagg 276*59599516SKenneth E. Jansen end if 277*59599516SKenneth E. Jansen enddo 278*59599516SKenneth E. Jansen 279*59599516SKenneth E. Jansen end subroutine !ramg_PCG 280*59599516SKenneth E. Jansen 281*59599516SKenneth E. Jansen!******************************************* 282*59599516SKenneth E. Jansen! ramg Interface 283*59599516SKenneth E. Jansen! Interface for libLES 284*59599516SKenneth E. Jansen!******************************************* 285*59599516SKenneth E. Jansen subroutine ramg_interface( 286*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,flowDiag, 287*59599516SKenneth E. Jansen & mcgR,mcgZ, 288*59599516SKenneth E. Jansen & ilwork,BC,iBC,iper 289*59599516SKenneth E. Jansen & ) 290*59599516SKenneth E. Jansen use ramg_data 291*59599516SKenneth E. Jansen include "common.h" 292*59599516SKenneth E. Jansen include "mpif.h" 293*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg+1) :: colm 294*59599516SKenneth E. Jansen integer,intent(in),dimension(nnz_tot) :: rowp 295*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK 296*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP 297*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg) :: mcgR 298*59599516SKenneth E. Jansen real(kind=8),intent(inout),dimension(nshg) :: mcgZ 299*59599516SKenneth E. Jansen integer, intent(in), dimension(nlwork) :: ilwork 300*59599516SKenneth E. Jansen integer, intent(in),dimension(nshg) :: iBC,iper 301*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC 302*59599516SKenneth E. Jansen real(kind=8),dimension(nshg) :: myR 303*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,4) :: flowDiag 304*59599516SKenneth E. Jansen 305*59599516SKenneth E. Jansen integer::i,j,k,amgmode 306*59599516SKenneth E. Jansen logical :: precflag 307*59599516SKenneth E. Jansen !character*10 :: fname 308*59599516SKenneth E. Jansen 309*59599516SKenneth E. Jansen mcgZ = 0 310*59599516SKenneth E. Jansen amgmode = 11 311*59599516SKenneth E. Jansen 312*59599516SKenneth E. Jansen ramg_winct = mod(ramg_winct+1,2) 313*59599516SKenneth E. Jansen 314*59599516SKenneth E. Jansen if (irun_amg_prec.eq.0) then 315*59599516SKenneth E. Jansen precflag = .false. 316*59599516SKenneth E. Jansen else if ((irun_amg_prec.ge.2).and.(ramg_flag.eq.1)) then 317*59599516SKenneth E. Jansen ! first solve attempt with CG if using smart solver 318*59599516SKenneth E. Jansen precflag=.false. 319*59599516SKenneth E. Jansen else 320*59599516SKenneth E. Jansen precflag = .true. 321*59599516SKenneth E. Jansen endif 322*59599516SKenneth E. Jansen 323*59599516SKenneth E. Jansen if (precflag) then !.and.(ramg_redo.ne.0)) then 324*59599516SKenneth E. Jansen myR = mcgR*amg_ppeDiag(1)%p 325*59599516SKenneth E. Jansen call ramg_driver(colm,rowp,lhsK,lhsP, 326*59599516SKenneth E. Jansen & nshg,nnz_tot,nflow, 327*59599516SKenneth E. Jansen & myR, 328*59599516SKenneth E. Jansen & nlwork,ilwork,ndofBC,BC,iBC,iper,mcgZ,1, 329*59599516SKenneth E. Jansen & ramg_eps,amgmode) 330*59599516SKenneth E. Jansen mcgZ = mcgZ*amg_ppeDiag(1)%p 331*59599516SKenneth E. Jansen else 332*59599516SKenneth E. Jansen mcgZ = mcgR 333*59599516SKenneth E. Jansen endif 334*59599516SKenneth E. Jansen 335*59599516SKenneth E. Jansen end subroutine ramg_interface 336*59599516SKenneth E. Jansen 337*59599516SKenneth E. Jansen!******************************************* 338*59599516SKenneth E. Jansen! ramg Driver 339*59599516SKenneth E. Jansen! 1. extract PPE / system 340*59599516SKenneth E. Jansen! 2. call coarsening 341*59599516SKenneth E. Jansen! 3. solve 342*59599516SKenneth E. Jansen!******************************************** 343*59599516SKenneth E. Jansen 344*59599516SKenneth E. Jansen !************************************* 345*59599516SKenneth E. Jansen ! ramg_driver 346*59599516SKenneth E. Jansen ! Input: global matrix,# of systems 347*59599516SKenneth E. Jansen ! control params 348*59599516SKenneth E. Jansen ! Output: solution 349*59599516SKenneth E. Jansen !************************************* 350*59599516SKenneth E. Jansen subroutine ramg_driver( 351*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP, 352*59599516SKenneth E. Jansen & nshg,nnz_tot,nflow, 353*59599516SKenneth E. Jansen & pperhs, 354*59599516SKenneth E. Jansen & nlwork,ilwork,ndofBC,BC,iBC,iper, 355*59599516SKenneth E. Jansen & ramg_sol,n_sol, 356*59599516SKenneth E. Jansen & amg_eps,amg_mode 357*59599516SKenneth E. Jansen & ) 358*59599516SKenneth E. Jansen 359*59599516SKenneth E. Jansen use ramg_data 360*59599516SKenneth E. Jansen implicit none 361*59599516SKenneth E. Jansen 362*59599516SKenneth E. Jansen !***********parameters************** 363*59599516SKenneth E. Jansen !the matrix 364*59599516SKenneth E. Jansen integer,intent(in) :: nshg 365*59599516SKenneth E. Jansen integer,intent(in) :: nnz_tot 366*59599516SKenneth E. Jansen integer,intent(in) :: nflow 367*59599516SKenneth E. Jansen !the matrix 368*59599516SKenneth E. Jansen integer,intent(in),dimension(nshg+1) :: colm 369*59599516SKenneth E. Jansen integer,intent(in),dimension(nnz_tot) :: rowp 370*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(9,nnz_tot) :: lhsK 371*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(4,nnz_tot) :: lhsP 372*59599516SKenneth E. Jansen !the forcing term, rhs 373*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg) :: pperhs 374*59599516SKenneth E. Jansen ! the boundary info 375*59599516SKenneth E. Jansen integer, intent(in) :: nlwork 376*59599516SKenneth E. Jansen integer, intent(in), dimension(nlwork) :: ilwork 377*59599516SKenneth E. Jansen integer, intent(in) :: ndofBC 378*59599516SKenneth E. Jansen integer, intent(in),dimension(nshg) :: iBC,iper 379*59599516SKenneth E. Jansen real(kind=8),intent(in),dimension(nshg,ndofBC) :: BC 380*59599516SKenneth E. Jansen !the output solution 381*59599516SKenneth E. Jansen integer, intent(in) :: n_sol 382*59599516SKenneth E. Jansen real(kind=8),intent(inout),dimension(nshg,n_sol) :: ramg_sol 383*59599516SKenneth E. Jansen 384*59599516SKenneth E. Jansen !the tolerance 385*59599516SKenneth E. Jansen real(kind=8),intent(in) :: amg_eps 386*59599516SKenneth E. Jansen !control run mode 387*59599516SKenneth E. Jansen integer,intent(inout) :: amg_mode 388*59599516SKenneth E. Jansen ! my AMG parameter 389*59599516SKenneth E. Jansen !*********parameters end************** 390*59599516SKenneth E. Jansen 391*59599516SKenneth E. Jansen !*****local variables***************** 392*59599516SKenneth E. Jansen real(kind=8) :: cpusec(10) 393*59599516SKenneth E. Jansen 394*59599516SKenneth E. Jansen call cpu_time(cpusec(1)) 395*59599516SKenneth E. Jansen 396*59599516SKenneth E. Jansen if (amg_mode .eq. 11 ) then ! solve PPE les Precondition ramg 397*59599516SKenneth E. Jansen call cpu_time(cpusec(2)) 398*59599516SKenneth E. Jansen call ramg_V_cycle(ramg_sol,pperhs,1, 399*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 400*59599516SKenneth E. Jansen call ramg_G_cycle(ramg_sol,pperhs,1, 401*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 402*59599516SKenneth E. Jansen call cpu_time(cpusec(3)) 403*59599516SKenneth E. Jansen ramg_time(4)=ramg_time(4)+cpusec(3)-cpusec(2) 404*59599516SKenneth E. Jansen else if (amg_mode .eq. 3) then ! solve PPE CG 405*59599516SKenneth E. Jansen call cpu_time(cpusec(2)) 406*59599516SKenneth E. Jansen call ramg_PCG(ramg_sol,pperhs,1, 407*59599516SKenneth E. Jansen & colm,rowp,lhsK,lhsP,ilwork,BC,iBC,iper) 408*59599516SKenneth E. Jansen call cpu_time(cpusec(3)) 409*59599516SKenneth E. Jansen ramg_time(4)=ramg_time(4)+cpusec(3)-cpusec(2) 410*59599516SKenneth E. Jansen end if 411*59599516SKenneth E. Jansen 412*59599516SKenneth E. Jansen return 413*59599516SKenneth E. Jansen 414*59599516SKenneth E. Jansen end subroutine ramg_driver 415*59599516SKenneth E. Jansen!******************************************* 416*59599516SKenneth E. Jansen! <EOF> ramg Driver 417*59599516SKenneth E. Jansen!******************************************* 418