1C Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2C the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3C reserved. See files LICENSE and NOTICE for details. 4C 5C This file is part of CEED, a collection of benchmarks, miniapps, software 6C libraries and APIs for efficient high-order finite element and spectral 7C element discretizations for exascale applications. For more information and 8C source code availability see http://github.com/ceed. 9C 10C The CEED research is supported by the Exascale Computing Project (17-SC-20-SC) 11C a collaborative effort of two U.S. Department of Energy organizations (Office 12C of Science and the National Nuclear Security Administration) responsible for 13C the planning and preparation of a capable exascale ecosystem, including 14C software, applications, hardware, advanced system engineering and early 15C testbed platforms, in support of the nation's exascale computing imperative. 16 17C> @file 18C> Mass and diffusion operators examples using Nek5000 19C TESTARGS -c {ceed_resource} -e bp1 -n 1 -b 4 -test 20 21C----------------------------------------------------------------------- 22 subroutine masssetupf(ctx,q,u1,u2,u3,u4,u5,u6,u7, 23 $ u8,u9,u10,u11,u12,u13,u14,u15,u16,v1,v2,v3,v4,v5,v6,v7,v8, 24 $ v9,v10,v11,v12,v13,v14,v15,v16,ierr) 25C Set up mass operator 26C Input: u1,u2,u3,q Output: v1,v2,ierr 27 integer q,ierr 28 real*8 ctx(1) 29 real*8 u1(3*q) 30 real*8 u2(9*q) 31 real*8 u3(q) 32 real*8 v1(q) 33 real*8 v2(q) 34 real*8 a11,a12,a13,a21,a22,a23,a31,a32,a33 35 real*8 g11,g12,g13,g21,g22,g23,g31,g32,g33 36 real*8 jacmq 37 38C Quadrature Point Loop 39 do i=1,q 40 a11=u2(i+q*0) 41 a12=u2(i+q*3) 42 a13=u2(i+q*6) 43 44 a21=u2(i+q*1) 45 a22=u2(i+q*4) 46 a23=u2(i+q*7) 47 48 a31=u2(i+q*2) 49 a32=u2(i+q*5) 50 a33=u2(i+q*8) 51 52 g11 = (a22*a33-a23*a32) 53 g12 = (a13*a32-a33*a12) 54 g13 = (a12*a23-a22*a13) 55 56 g21 = (a23*a31-a21*a33) 57 g22 = (a11*a33-a31*a13) 58 g23 = (a13*a21-a23*a11) 59 60 g31 = (a21*a32-a22*a31) 61 g32 = (a12*a31-a32*a11) 62 g33 = (a11*a22-a21*a12) 63 64 jacmq = a11*g11+a21*g12+a31*g13 65 66C Rho 67 v1(i)=u3(i)*jacmq 68 69C RHS 70 v2(i)=u3(i)*jacmq 71 $ *dsqrt(u1(i+q*0)*u1(i+q*0) 72 $ +u1(i+q*1)*u1(i+q*1) 73 $ +u1(i+q*2)*u1(i+q*2)) 74 enddo 75 76 ierr=0 77 end 78C----------------------------------------------------------------------- 79 subroutine massf(ctx,q,u1,u2,u3,u4,u5,u6,u7, 80 $ u8,u9,u10,u11,u12,u13,u14,u15,u16,v1,v2,v3,v4,v5,v6,v7,v8, 81 $ v9,v10,v11,v12,v13,v14,v15,v16,ierr) 82C Apply mass operator 83C Input: u1,u2,q Output: v1,ierr 84 integer q,ierr 85 real*8 ctx(1) 86 real*8 u1(q) 87 real*8 u2(q) 88 real*8 v1(q) 89 90C Quadrature Point Loop 91 do i=1,q 92 v1(i)=u2(i)*u1(i) 93 enddo 94 95 ierr=0 96 end 97C----------------------------------------------------------------------- 98 subroutine diffsetupf(ctx,q,u1,u2,u3,u4,u5,u6,u7, 99 $ u8,u9,u10,u11,u12,u13,u14,u15,u16,v1,v2,v3,v4,v5,v6,v7,v8, 100 $ v9,v10,v11,v12,v13,v14,v15,v16,ierr) 101C Set up diffusion operator 102C Input: u1,u2,u3,q Output: v1,v2,ierr 103 integer q,ierr 104 real*8 ctx(1) 105 real*8 u1(3*q) 106 real*8 u2(9*q) 107 real*8 u3(q) 108 real*8 v1(6*q) 109 real*8 v2(q) 110 real*8 a11,a12,a13,a21,a22,a23,a31,a32,a33 111 real*8 g11,g12,g13,g21,g22,g23,g31,g32,g33 112 real*8 jacmq,scl 113 real*8 c(3),k(3) 114 115C Quadrature Point Loop 116 do i=1,q 117 pi = 3.14159265358979323846 118 119 c(1)=0. 120 c(2)=1. 121 c(3)=2. 122 k(1)=1. 123 k(2)=2. 124 k(3)=3. 125 126 a11=u2(i+q*0) 127 a12=u2(i+q*3) 128 a13=u2(i+q*6) 129 130 a21=u2(i+q*1) 131 a22=u2(i+q*4) 132 a23=u2(i+q*7) 133 134 a31=u2(i+q*2) 135 a32=u2(i+q*5) 136 a33=u2(i+q*8) 137 138 g11 = (a22*a33-a23*a32) 139 g12 = (a13*a32-a33*a12) 140 g13 = (a12*a23-a22*a13) 141 142 g21 = (a23*a31-a21*a33) 143 g22 = (a11*a33-a31*a13) 144 g23 = (a13*a21-a23*a11) 145 146 g31 = (a21*a32-a22*a31) 147 g32 = (a12*a31-a32*a11) 148 g33 = (a11*a22-a21*a12) 149 150 jacmq = a11*g11+a21*g12+a31*g13 151 152 scl = u3(i)/jacmq 153 154C Geometric factors 155C Stored in Voigt convention 156C 0 5 4 157C 5 1 3 158C 4 3 2 159 v1(i+0*q) = scl*(g11*g11+g12*g12+g13*g13) ! Grr 160 v1(i+1*q) = scl*(g21*g21+g22*g22+g23*g23) ! Gss 161 v1(i+2*q) = scl*(g31*g31+g32*g32+g33*g33) ! Gtt 162 v1(i+3*q) = scl*(g21*g31+g22*g32+g23*g33) ! Gst 163 v1(i+4*q) = scl*(g11*g31+g12*g32+g13*g33) ! Grt 164 v1(i+5*q) = scl*(g11*g21+g12*g22+g13*g23) ! Grs 165 166C RHS 167 v2(i) = u3(i)*jacmq*pi*pi 168 $ *dsin(pi*(c(1)+k(1)*u1(i+0*q))) 169 $ *dsin(pi*(c(2)+k(2)*u1(i+1*q))) 170 $ *dsin(pi*(c(3)+k(3)*u1(i+2*q))) 171 $ *(k(1)*k(1)+k(2)*k(2)+k(3)*k(3)) 172 173 enddo 174 175 ierr=0 176 end 177C----------------------------------------------------------------------- 178 subroutine diffusionf(ctx,q,u1,u2,u3,u4,u5,u6,u7, 179 $ u8,u9,u10,u11,u12,u13,u14,u15,u16,v1,v2,v3,v4,v5,v6,v7,v8, 180 $ v9,v10,v11,v12,v13,v14,v15,v16,ierr) 181C Apply diffusion operator 182C Input: u1,u2,q Output: v1,ierr 183 integer q,ierr 184 real*8 ctx(1) 185 real*8 u1(3*q) 186 real*8 u2(6*q) 187 real*8 v1(3*q) 188 189C Quadrature Point Loop 190 do i=1,q 191 v1(i+0*q)= 192 $ u2(i+0*q)*u1(i)+u2(i+5*q)*u1(i+q)+u2(i+4*q)*u1(i+2*q) 193 v1(i+1*q)= 194 $ u2(i+5*q)*u1(i)+u2(i+1*q)*u1(i+q)+u2(i+3*q)*u1(i+2*q) 195 v1(i+2*q)= 196 $ u2(i+4*q)*u1(i)+u2(i+3*q)*u1(i+q)+u2(i+2*q)*u1(i+2*q) 197 enddo 198 199 ierr=0 200 end 201C----------------------------------------------------------------------- 202 subroutine set_h2_as_rhoJac_GL(h2,bmq,nxq) 203C Set h2 as rhoJac 204C Input: bmq,nxq Output: h2 205 include 'SIZE' 206 real*8 h2(1),bmq(1) 207 208 common /ctmp77/ wd(lxd),zd(lxd) 209 integer e,i,L 210 211 call zwgl(zd,wd,nxq) ! nxq = number of points 212 213 q = 1.0 ! Later, this can be a function of position... 214 215 L = 0 216 do e=1,lelt 217 do i=1,nxq**ldim 218 L=L+1 219 h2(L) = q*bmq(L) 220 enddo 221 enddo 222 223 return 224 end 225C----------------------------------------------------------------------- 226 subroutine dist_fld_h1(e) 227C Set distance initial condition for BP1 228C Input: Output: e 229 include 'SIZE' 230 include 'TOTAL' 231 real*8 x,y,z 232 real*8 e(1) 233 234 n=lx1*ly1*lz1*nelt 235 236 do i=1,n 237 x=xm1(i,1,1,1) 238 y=ym1(i,1,1,1) 239 z=zm1(i,1,1,1) 240 241 e(i) = dsqrt(x*x+y*y+z*z) 242 enddo 243 244 call dsavg(e) ! This is requisite for random fields 245 246 return 247 end 248C----------------------------------------------------------------------- 249 subroutine sin_fld_h1(e) 250C Set sine initial condition for BP3 251C Input: Output: e 252 include 'SIZE' 253 include 'TOTAL' 254 real*8 x,y,z 255 real*8 e(1) 256 real*8 c(3),k(3) 257 258 n=lx1*ly1*lz1*nelt 259 pi = 3.14159265358979323846 260 261 c(1)=0. 262 c(2)=1. 263 c(3)=2. 264 k(1)=1. 265 k(2)=2. 266 k(3)=3. 267 268 do i=1,n 269 x=xm1(i,1,1,1) 270 y=ym1(i,1,1,1) 271 z=zm1(i,1,1,1) 272 273 e(i) = dsin(pi*(c(1)+k(1)*x)) 274 & *dsin(pi*(c(2)+k(2)*y)) 275 & *dsin(pi*(c(3)+k(3)*z)) 276 277 enddo 278 279 call dsavg(e) ! This is requisite for random fields 280 281 return 282 end 283C----------------------------------------------------------------------- 284 subroutine uservp(ix,iy,iz,eg) ! set variable properties 285 include 'SIZE' 286 include 'TOTAL' 287 include 'NEKUSE' 288 integer e,f,eg 289C e = gllel(eg) 290 291 udiff = 0.0 292 utrans = 0.0 293 294 return 295 end 296C----------------------------------------------------------------------- 297 subroutine userf(ix,iy,iz,eg) ! set acceleration term 298C 299C Note: this is an acceleration term, NOT a force! 300C Thus, ffx will subsequently be multiplied by rho(x,t). 301C 302 include 'SIZE' 303 include 'TOTAL' 304 include 'NEKUSE' 305 integer e,f,eg 306C e = gllel(eg) 307 308 ffx = 0.0 309 ffy = 0.0 310 ffz = 0.0 311 312 return 313 end 314C----------------------------------------------------------------------- 315 subroutine userq(i,j,k,eg) ! set source term 316 include 'SIZE' 317 include 'TOTAL' 318 include 'NEKUSE' 319 integer e,f,eg 320 e = gllel(eg) 321 322 qvol = 0 323 324 return 325 end 326C----------------------------------------------------------------------- 327 subroutine userbc(ix,iy,iz,f,eg) ! set up boundary conditions 328C NOTE ::: This subroutine MAY NOT be called by every process 329 include 'SIZE' 330 include 'TOTAL' 331 include 'NEKUSE' 332 integer e,f,eg 333 334 ux = 0.0 335 uy = 0.0 336 uz = 0.0 337 temp = 0.0 338 339 return 340 end 341C----------------------------------------------------------------------- 342 subroutine useric(ix,iy,iz,eg) ! set up initial conditions 343 include 'SIZE' 344 include 'TOTAL' 345 include 'NEKUSE' 346 integer e,f,eg 347 348 ux = 0.0 349 uy = 0.0 350 uz = 0.0 351 temp = 0.0 352 353 return 354 end 355C----------------------------------------------------------------------- 356 subroutine usrdat ! This routine to modify element vertices 357 include 'SIZE' 358 include 'TOTAL' 359 360 return 361 end 362C----------------------------------------------------------------------- 363 subroutine usrdat2 ! This routine to modify mesh coordinates 364 include 'SIZE' 365 include 'TOTAL' 366 367 x0 = 0 368 x1 = 1 369 call rescale_x(xm1,x0,x1) 370 call rescale_x(ym1,x0,x1) 371 call rescale_x(zm1,x0,x1) 372 373 param(59)=1 ! Force Nek to use the "deformed element" formulation 374 375 return 376 end 377C----------------------------------------------------------------------- 378 subroutine usrdat3 379 include 'SIZE' 380 include 'TOTAL' 381 382 return 383 end 384C----------------------------------------------------------------------- 385 subroutine xmask1 (r1,h2,nel) 386C Apply zero Dirichlet boundary conditions 387C Input: h2,nel Output: r1 388 include 'SIZE' 389 include 'TOTAL' 390 real*8 r1(1),h2(1) 391 392 n=nx1*ny1*nz1*nel 393 do i=1,n 394 r1(i)=r1(i)*h2(i) 395 enddo 396 397 return 398 end 399C----------------------------------------------------------------------- 400 function glrdif(x,y,n) 401C Compute Linfty norm of (x-y) 402C Input: x,y Output: n 403 real*8 x(n),y(n) 404 405 dmx=0 406 xmx=0 407 ymx=0 408 409 do i=1,n 410 diff=abs(x(i)-y(i)) 411 dmx =max(dmx,diff) 412 xmx =max(xmx,x(i)) 413 ymx =max(ymx,y(i)) 414 enddo 415 416 xmx = max(xmx,ymx) 417 dmx = glmax(dmx,1) ! max across processors 418 xmx = glmax(xmx,1) 419 420 if (xmx.gt.0) then 421 glrdif = dmx/xmx 422 else 423 glrdif = -dmx ! Negative indicates something strange 424 endif 425 426 return 427 end 428C----------------------------------------------------------------------- 429 subroutine loc_grad3(ur,us,ut,u,N,D,Dt) 430C 3D transpose of local gradient 431C Input: u,N,D,Dt Output: ur,us,ut 432 real*8 ur(0:N,0:N,0:N),us(0:N,0:N,0:N),ut(0:N,0:N,0:N) 433 real*8 u (0:N,0:N,0:N) 434 real*8 D (0:N,0:N),Dt(0:N,0:N) 435 436 m1 = N+1 437 m2 = m1*m1 438 439 call mxm(D ,m1,u(0,0,0),m1,ur,m2) 440 do k=0,N 441 call mxm(u(0,0,k),m1,Dt,m1,us(0,0,k),m1) 442 enddo 443 call mxm(u(0,0,0),m2,Dt,m1,ut,m1) 444 445 return 446 end 447c----------------------------------------------------------------------- 448 subroutine loc_grad3t(u,ur,us,ut,N,D,Dt,w) 449C 3D transpose of local gradient 450C Input: ur,us,ut,N,D,Dt Output: u 451 real*8 u (0:N,0:N,0:N) 452 real*8 ur(0:N,0:N,0:N),us(0:N,0:N,0:N),ut(0:N,0:N,0:N) 453 real*8 D (0:N,0:N),Dt(0:N,0:N) 454 real*8 w (0:N,0:N,0:N) 455 456 m1 = N+1 457 m2 = m1*m1 458 m3 = m1*m1*m1 459 460 call mxm(Dt,m1,ur,m1,u(0,0,0),m2) 461 do k=0,N 462 call mxm(us(0,0,k),m1,D ,m1,w(0,0,k),m1) 463 enddo 464 call add2(u(0,0,0),w,m3) 465 call mxm(ut,m2,D ,m1,w,m1) 466 call add2(u(0,0,0),w,m3) 467 468 return 469 end 470C----------------------------------------------------------------------- 471 subroutine geodatq(gf,bmq,w3mq,nxq) 472C Routine to generate elemental geometric matrices on mesh 1 473C (Gauss-Legendre Lobatto mesh). 474 include 'SIZE' 475 include 'TOTAL' 476 477 parameter (lg=3+3*(ldim-2),lzq=lx1+1,lxyd=lzq**ldim) 478 479 real*8 gf(lg,nxq**ldim,lelt),bmq(nxq**ldim,lelt),w3mq(nxq,nxq,nxq) 480 481 common /ctmp1/ xr(lxyd),xs(lxyd),xt(lxyd) 482 common /sxrns/ yr(lxyd),ys(lxyd),yt(lxyd) 483 $ , zr(lxyd),zs(lxyd),zt(lxyd) 484 485 common /ctmp77/ wd(lxd),zd(lxd) 486 common /dxmfine/ dxmq(lzq,lzq),dxtmq(lzq,lzq) 487 488 integer e 489 real*8 tmp(lxyd) 490 real*8 a11,a12,a13,a21,a22,a23,a31,a32,a33 491 real*8 g11,g12,g13,g21,g22,g23,g31,g32,g33 492 real*8 jacmq 493 494 if (nxq.gt.lzq) call exitti('ABORT: recompile with lzq=$',nxq) 495 496 call zwgl (zd,wd,lzq) ! nxq = number of points 497 call gen_dgl (dxmq,dxtmq,lzq,lzq,tmp) 498 499 do k=1,nxq 500 do j=1,nxq 501 do i=1,nxq 502 w3mq(i,j,k) = wd(i)*wd(j)*wd(k) 503 enddo 504 enddo 505 enddo 506 507 nxyzq = nxq**ldim 508 nxqm1 = lzq-1 509 510 do e=1,nelt 511 call intp_rstd (tmp,xm1(1,1,1,e),lx1,lzq,if3d,0) ! 0-->Fwd interpolation 512 call loc_grad3 (xr,xs,xt,tmp,nxqm1,dxmq,dxtmq) 513 514 call intp_rstd (tmp,ym1(1,1,1,e),lx1,lzq,if3d,0) 515 call loc_grad3 (yr,ys,yt,tmp,nxqm1,dxmq,dxtmq) 516 517 call intp_rstd (tmp,zm1(1,1,1,e),lx1,lzq,if3d,0) 518 call loc_grad3 (zr,zs,zt,tmp,nxqm1,dxmq,dxtmq) 519 520 do i=1,nxyzq 521 a11 = xr(i) 522 a12 = xs(i) 523 a13 = xt(i) 524 525 a21 = yr(i) 526 a22 = ys(i) 527 a23 = yt(i) 528 529 a31 = zr(i) 530 a32 = zs(i) 531 a33 = zt(i) 532 533 g11 = (a22*a33-a23*a32) 534 g12 = (a13*a32-a33*a12) 535 g13 = (a12*a23-a22*a13) 536 537 g21 = (a23*a31-a21*a33) 538 g22 = (a11*a33-a31*a13) 539 g23 = (a13*a21-a23*a11) 540 541 g31 = (a21*a32-a22*a31) 542 g32 = (a12*a31-a32*a11) 543 g33 = (a11*a22-a21*a12) 544 545 jacmq = a11*g11+a21*g12+a31*g13 546 547 bmq(i,e) = w3mq(i,1,1)*jacmq 548 scale = w3mq(i,1,1)/jacmq 549 550 gf(1,i,e) = scale*(g11*g11+g12*g12+g13*g13) ! Grr 551 gf(2,i,e) = scale*(g11*g21+g12*g22+g13*g23) ! Grs 552 gf(3,i,e) = scale*(g11*g31+g12*g32+g13*g33) ! Grt 553 gf(4,i,e) = scale*(g21*g21+g22*g22+g23*g23) ! Gss 554 gf(5,i,e) = scale*(g21*g31+g22*g32+g23*g33) ! Gst 555 gf(6,i,e) = scale*(g31*g31+g32*g32+g33*g33) ! Gtt 556 557 enddo 558 enddo 559 560 return 561 end 562C----------------------------------------------------------------------- 563 subroutine setprecn_bp1 (d,h1,h2) 564C Generate diagonal preconditioner for Helmholtz operator 565C Input: h1,h2 Output: d 566 include 'SIZE' 567 include 'TOTAL' 568 569 parameter (lxyz=lx1*ly1*lz1,lg=3+3*(ldim-2)) 570 571 real*8 d(lx1,ly1,lz1,lelt),h1(lxyz,lelt),h2(lxyz,lelt) 572 integer e 573 574 real*8 gf(lg,lx1,ly1,lz1,lelt) ! Equivalence new gf() data 575 equivalence (gf,g1m1) ! layout to g1m1...g6m1 576 577 real*8 ysm1(ly1) 578 579 nel = nelfld(ifield) 580 n = nel*lx1*ly1*lz1 581 nxyz = lx1*ly1*lz1 582 583 call copy (d,bm1,n) ! Mass matrix preconditioning full mass matrix 584 call dssum (d,nx1,ny1,nz1) 585 call invcol1 (d,n) 586 return 587 588 call dsset(lx1,ly1,lz1) 589 590 do 1000 e=1,nel 591 592 call rzero(d(1,1,1,e),nxyz) 593 594 do 320 iz=1,lz1 595 do 320 iy=1,ly1 596 do 320 ix=1,lx1 597 do 320 iq=1,lx1 598 d(ix,iy,iz,e) = d(ix,iy,iz,e) 599 $ + gf(1,iq,iy,iz,e) * dxm1(iq,ix)**2 600 $ + gf(2,ix,iq,iz,e) * dxm1(iq,iy)**2 601 $ + gf(3,ix,iy,iq,e) * dxm1(iq,iz)**2 602 320 continue 603C 604C Add cross terms if element is deformed. 605C 606 if (lxyz.gt.0) then 607 608 do i2=1,ly1,ly1-1 609 do i1=1,lx1,lx1-1 610 d(1,i1,i2,e) = d(1,i1,i2,e) 611 $ + gf(4,1,i1,i2,e) * dxtm1(1,1)*dytm1(i1,i1) 612 $ + gf(5,1,i1,i2,e) * dxtm1(1,1)*dztm1(i2,i2) 613 d(lx1,i1,i2,e) = d(lx1,i1,i2,e) 614 $ + gf(4,lx1,i1,i2,e) * dxtm1(lx1,lx1)*dytm1(i1,i1) 615 $ + gf(5,lx1,i1,i2,e) * dxtm1(lx1,lx1)*dztm1(i2,i2) 616 d(i1,1,i2,e) = d(i1,1,i2,e) 617 $ + gf(4,i1,1,i2,e) * dytm1(1,1)*dxtm1(i1,i1) 618 $ + gf(6,i1,1,i2,e) * dytm1(1,1)*dztm1(i2,i2) 619 d(i1,ly1,i2,e) = d(i1,ly1,i2,e) 620 $ + gf(4,i1,ly1,i2,e) * dytm1(ly1,ly1)*dxtm1(i1,i1) 621 $ + gf(6,i1,ly1,i2,e) * dytm1(ly1,ly1)*dztm1(i2,i2) 622 d(i1,i2,1,e) = d(i1,i2,1,e) 623 $ + gf(5,i1,i2,1,e) * dztm1(1,1)*dxtm1(i1,i1) 624 $ + gf(6,i1,i2,1,e) * dztm1(1,1)*dytm1(i2,i2) 625 d(i1,i2,lz1,e) = d(i1,i2,lz1,e) 626 $ + gf(5,i1,i2,lz1,e) * dztm1(lz1,lz1)*dxtm1(i1,i1) 627 $ + gf(6,i1,i2,lz1,e) * dztm1(lz1,lz1)*dytm1(i2,i2) 628 629 enddo 630 enddo 631 endif 632 633 do i=1,lxyz 634 d(i,1,1,e)=d(i,1,1,e)*h1(i,e)+h2(i,e)*bm1(i,1,1,e) 635 enddo 636 637 1000 continue ! element loop 638 639C If axisymmetric, add a diagonal term in the radial direction (ISD=2) 640 641 if (ifaxis.and.(isd.eq.2)) then 642 do 1200 e=1,nel 643 if (ifrzer(e)) call mxm(ym1(1,1,1,e),lx1,datm1,ly1,ysm1,1) 644 k=0 645 do 1190 j=1,ly1 646 do 1190 i=1,lx1 647 k=k+1 648 if (ym1(i,j,1,e).ne.0.) then 649 term1 = bm1(i,j,1,e)/ym1(i,j,1,e)**2 650 if (ifrzer(e)) then 651 term2 = wxm1(i)*wam1(1)*dam1(1,j) 652 $ *jacm1(i,1,1,e)/ysm1(i) 653 else 654 term2 = 0. 655 endif 656 d(i,j,1,e) = d(i,j,1,e) + h1(k,e)*(term1+term2) 657 endif 658 1190 continue 659 1200 continue 660 661 endif 662 call dssum (d,nx1,ny1,nz1) 663 call invcol1 (d,n) 664 665 if (nio.eq.0) write(6,1) n,d(1,1,1,1),h1(1,1),h2(1,1),bm1(1,1,1,1) 666 1 format(i9,1p4e12.4,' diag prec') 667 668 return 669 end 670C----------------------------------------------------------------------- 671 subroutine setprecn_bp3 (d,h1,h2) 672C Generate dummy diagonal preconditioner for Helmholtz operator 673C Input: h1,h2 Output: d 674 include 'SIZE' 675 include 'TOTAL' 676 677 parameter (n=lx1*ly1*lz1*lelt) 678 real*8 d(n),h1(n),h2(n) 679 680 call rone (d,n) 681 682 return 683 end 684C----------------------------------------------------------------------- 685 subroutine userchk 686 include 'SIZE' 687 include 'TOTAL' 688 689 integer bp 690 691 call get_bp(bp) 692 693 if (bp==1) then 694 if (istep.gt.0) call bp1 695 elseif (bp==3) then 696 if (istep.gt.0) call bp3 697 else 698 write(6,*) "INVALID BP SPECIFICED" 699 endif 700 701 return 702 end 703C----------------------------------------------------------------------- 704 subroutine bp1 705C Solution to BP1 using libCEED 706 include 'SIZE' 707 include 'TOTAL' 708 include 'CTIMER' ! ifsync 709 include 'FDMH1' 710 include 'ceedf.h' 711 712 parameter (lzq=lx1+1) 713 parameter (lx=lx1*ly1*lz1,lg=3+3*(ldim-2),lq=lzq**ldim) 714 common /bpgfactors/ gf(lg*lq,lelt),bmq(lq,lelt),w3mq(lq) 715 716 parameter (lt=lx1*ly1*lz1*lelt) 717 parameter (ld=lxd*lyd*lzd*lelt) 718 common /vcrns/ u1(lt),r1(lt),r2(lt),r3(lt) 719 common /vcrny/ e1(lt) 720 common /vcrvh/ h1(lt),h2(lx*lelt),pap(3) 721 real*8 coords(ldim*lx*lelt) 722 723 logical ifh3 724 integer*8 nnode 725 integer ceed,err,test 726 character*64 spec 727 728 integer p,q,ncompx,ncompu,enode,lnode 729 integer vec_p1,vec_ap1,vec_qdata,vec_coords,vec_rhs 730 integer erstrctu,erstrctx,erstrctw 731 integer basisu,basisx 732 integer qf_mass,qf_setup 733 integer op_mass,op_setup 734 real*8 x,y,z 735 integer*8 offset 736 737 external massf,masssetupf 738 739 ifield = 1 740 nxq = nx1+1 741 n = nx1*ny1*nz1*nelt 742 743 ifsync = .false. 744 745C Set up coordinates 746 ii=0 747 do j=0,nelt-1 748 do i=1,lx 749 ii=ii+1 750 x = xm1(ii,1,1,1) 751 y = ym1(ii,1,1,1) 752 z = zm1(ii,1,1,1) 753 coords(i+0*lx+3*j*lx)=x 754 coords(i+1*lx+3*j*lx)=y 755 coords(i+2*lx+3*j*lx)=z 756 enddo 757 enddo 758 759C Init ceed library 760 call get_spec(spec) 761 call ceedinit(trim(spec)//char(0),ceed,err) 762 763 call get_test(test) 764 765C Set up Nek geometry data 766 call geodatq (gf,bmq,w3mq,nxq) ! Set up gf() arrays 767 call set_h2_as_rhoJac_GL(h2,bmq,nxq) 768 769C Set up true soln 770 call dist_fld_h1 (e1) 771 call copy (h1,e1,n) ! Save exact soln in h1 772 773C Set up solver parameters 774 tol = 1e-10 775 param(22) = tol 776 maxit = 100 777 778 call nekgsync() 779 780C Create ceed basis for mesh and computation 781 p=nx1 782 q=p+1 783 ncompu=1 784 ncompx=ldim 785 call ceedbasiscreatetensorh1lagrange(ceed,ndim,ndim,p,q, 786 $ ceed_gauss,basisx,err) 787 call ceedbasiscreatetensorh1lagrange(ceed,ndim,ncompu,p,q, 788 $ ceed_gauss,basisu,err) 789 790C Create ceed element restrictions for mesh and computation 791 enode=p**ldim 792 lnode=enode*nelt*ncompu 793 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 794 $ enode,lnode,ldim,erstrctx,err) 795 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 796 $ enode,lnode,ncompu,erstrctu,err) 797 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 798 $ q**ldim,nelt*q**ldim,1,erstrctw,err) 799 800C Create ceed vectors 801 call ceedvectorcreate(ceed,lnode,vec_p1,err) 802 call ceedvectorcreate(ceed,lnode,vec_ap1,err) 803 call ceedvectorcreate(ceed,lnode,vec_rhs,err) 804 call ceedvectorcreate(ceed,ldim*lx*nelt,vec_coords,err) 805 call ceedvectorcreate(ceed,nelt*q**ldim,vec_qdata,err) 806 807 offset=0 808 call ceedvectorsetarray(vec_coords,ceed_mem_host, 809 $ ceed_use_pointer,coords,offset,err) 810 811C Create ceed qfunctions for masssetupf and massf 812 call ceedqfunctioncreateinterior(ceed,1,masssetupf, 813 $ EXAMPLE_DIR 814 $ //'bps/bps.h:masssetupf',qf_setup,err) 815 call ceedqfunctionaddinput(qf_setup,'x',ncompx, 816 $ ceed_eval_interp,err) 817 call ceedqfunctionaddinput(qf_setup,'dx',ncompx*ldim, 818 $ ceed_eval_grad,err) 819 call ceedqfunctionaddinput(qf_setup,'weight',ncompu, 820 $ ceed_eval_weight,err) 821 call ceedqfunctionaddoutput(qf_setup,'qdata',ncompu, 822 $ ceed_eval_none,err) 823 call ceedqfunctionaddoutput(qf_setup,'rhs',ncompu, 824 $ ceed_eval_interp,err) 825 826 call ceedqfunctioncreateinterior(ceed,1,massf, 827 $ EXAMPLE_DIR 828 $ //'bps/bps.h:massf',qf_mass,err) 829 call ceedqfunctionaddinput(qf_mass,'u',ncompu, 830 $ ceed_eval_interp,err) 831 call ceedqfunctionaddinput(qf_mass,'qdata',ncompu, 832 $ ceed_eval_none,err) 833 call ceedqfunctionaddoutput(qf_mass,'v',ncompu, 834 $ ceed_eval_interp,err) 835 836C Create ceed operators 837 call ceedoperatorcreate(ceed,qf_setup, 838 $ ceed_qfunction_none,ceed_qfunction_none,op_setup,err) 839 call ceedoperatorsetfield(op_setup,'x',erstrctx, 840 $ basisx,ceed_vector_active,err) 841 call ceedoperatorsetfield(op_setup,'dx',erstrctx, 842 $ basisx,ceed_vector_active,err) 843 call ceedoperatorsetfield(op_setup,'weight',erstrctx, 844 $ basisx,ceed_vector_none,err) 845 call ceedoperatorsetfield(op_setup,'qdata',erstrctw, 846 $ ceed_basis_collocated,ceed_vector_active,err) 847 call ceedoperatorsetfield(op_setup,'rhs',erstrctu, 848 $ basisu,vec_rhs,err) 849 850 call ceedoperatorcreate(ceed,qf_mass, 851 $ ceed_qfunction_none,ceed_qfunction_none,op_mass,err) 852 call ceedoperatorsetfield(op_mass,'u',erstrctu, 853 $ basisu,ceed_vector_active,err) 854 call ceedoperatorsetfield(op_mass,'qdata',erstrctw, 855 $ ceed_basis_collocated,vec_qdata,err) 856 call ceedoperatorsetfield(op_mass,'v',erstrctu, 857 $ basisu,ceed_vector_active,err) 858 859C Compute setup data 860 call ceedvectorsetarray(vec_rhs,ceed_mem_host, 861 $ ceed_use_pointer,r1,offset,err) 862 call ceedoperatorapply(op_setup,vec_coords,vec_qdata, 863 $ ceed_request_immediate,err) 864 865C Set up true RHS 866 call dssum (r1,nx1,ny1,nz1) ! r1 867 868C Set up algebraic RHS with libCEED 869 call ceedvectorsetarray(vec_p1,ceed_mem_host, 870 $ ceed_use_pointer,h1,offset,err) 871 call ceedvectorsetarray(vec_ap1,ceed_mem_host, 872 $ ceed_use_pointer,r2,offset,err) 873 call ceedoperatorapply(op_mass,vec_p1,vec_ap1, 874 $ ceed_request_immediate,err) ! r2 = A_ceed*h1 875 call dssum (r2,nx1,ny1,nz1) 876 877C Set up algebraic RHS with Nek5000 878 call axhm1 (pap,r3,h1,h1,h2,'bp1') ! r3 = A_nek5k*h1 879 call dssum (r3,nx1,ny1,nz1) 880 881 call nekgsync() 882 883C Solve true RHS 884 if (nid.eq.0) write (6,*) "libCEED true RHS" 885 tstart = dnekclock() 886 call cggos(u1,r1,h1,h2,vmult,binvm1,tol,ceed,op_mass, 887 $ vec_p1,vec_ap1,maxit,'bp1') 888 tstop = dnekclock() 889 890C Output 891 telaps = (tstop-tstart) 892 maxits = maxit 893 er1 = glrdif(u1,e1,n) 894 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 895 896 if (test.eq.1.and.nid.eq.0) then 897 if (maxit>=100) then 898 write(6,*) "UNCONVERGED CG" 899 endif 900 if (dabs(er1)>5e-3) then 901 write(6,*) "ERROR IS TOO LARGE" 902 endif 903 endif 904 905 nx = nx1-1 906 nnode = nelgt ! nnodes 907 nnode = nnode*(nx**ldim) ! nnodes 908 nppp = nnode/np ! nnodes/proc 909 910 dofps = nnode/telaps ! DOF/sec - scalar form 911 titers = telaps/maxits ! time per iteration 912 tppp_s = titers/nppp ! time per iteraton per local point 913 914 if (nid.eq.0) write(6,1) 'case scalar:' 915 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,dofps,titers,tppp_s 916 917C Solve libCEED algebraic RHS 918 if (nid.eq.0) write (6,*) "libCEED algebraic RHS" 919 maxit = 100 920 tstart = dnekclock() 921 call cggos(u1,r2,h1,h2,vmult,binvm1,tol,ceed,op_mass, 922 $ vec_p1,vec_ap1,maxit,'bp1') 923 tstop = dnekclock() 924 925C Output 926 telaps = (tstop-tstart) 927 maxits = maxit 928 er1 = glrdif(u1,e1,n) 929 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 930 931 if (test.eq.1.and.nid.eq.0) then 932 if (maxit>=100) then 933 write(6,*) "UNCONVERGED CG" 934 endif 935 if (dabs(er1)>1e-5) then 936 write(6,*) "ERROR IS TOO LARGE" 937 endif 938 endif 939 940 nx = nx1-1 941 nnode = nelgt ! nnodes 942 nnode = nnode*(nx**ldim) ! nnodes 943 nppp = nnode/np ! nnodes/proc 944 945 dofps = nnode/telaps ! DOF/sec - scalar form 946 titers = telaps/maxits ! time per iteration 947 tppp_s = titers/nppp ! time per iteraton per local point 948 949 if (nid.eq.0) write(6,1) 'case scalar:' 950 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,dofps,titers,tppp_s 951 952C Solve Nek5000 algebraic RHS 953 if (nid.eq.0) write (6,*) "Nek5000 algebraic RHS" 954 maxit = 100 955 tstart = dnekclock() 956 call cggos(u1,r3,h1,h2,vmult,binvm1,tol,ceed,op_mass, 957 $ vec_p1,vec_ap1,maxit,'bp1') 958 tstop = dnekclock() 959 960C Output 961 telaps = (tstop-tstart) 962 maxits = maxit 963 er1 = glrdif(u1,e1,n) 964 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 965 966 if (test.eq.1.and.nid.eq.0) then 967 if (maxit>=100) then 968 write(6,*) "UNCONVERGED CG" 969 endif 970 if (dabs(er1)>1e-5) then 971 write(6,*) "ERROR IS TOO LARGE" 972 endif 973 endif 974 975 nx = nx1-1 976 nnode = nelgt ! nnodes 977 nnode = nnode*(nx**ldim) ! nnodes 978 nppp = nnode/np ! nnodes/proc 979 980 dofps = nnode/telaps ! DOF/sec - scalar form 981 titers = telaps/maxits ! time per iteration 982 tppp_s = titers/nppp ! time per iteraton per local point 983 984 if (nid.eq.0) write(6,1) 'case scalar:' 985 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,dofps,titers,tppp_s 986 987 1 format(a12,i7,i3,i7,i10,i14,i10,i4,1p4e13.5) 988 3 format(i3,i9,e12.4,1x,a8,i9) 989 990C Destroy ceed handles 991 call ceedvectordestroy(vec_p1,err) 992 call ceedvectordestroy(vec_ap1,err) 993 call ceedvectordestroy(vec_rhs,err) 994 call ceedvectordestroy(vec_qdata,err) 995 call ceedvectordestroy(vec_coords,err) 996 call ceedelemrestrictiondestroy(erstrctu,err) 997 call ceedelemrestrictiondestroy(erstrctx,err) 998 call ceedelemrestrictiondestroy(erstrctw,err) 999 call ceedbasisdestroy(basisu,err) 1000 call ceedbasisdestroy(basisx,err) 1001 call ceedqfunctiondestroy(qf_setup,err) 1002 call ceedqfunctiondestroy(qf_mass,err) 1003 call ceedoperatordestroy(op_setup,err) 1004 call ceedoperatordestroy(op_mass,err) 1005 call ceeddestroy(ceed,err) 1006 1007 return 1008 end 1009C----------------------------------------------------------------------- 1010 subroutine bp3 1011C Solution to BP3 using libCEED 1012 include 'SIZE' 1013 include 'TOTAL' 1014 include 'CTIMER' ! ifsync 1015 include 'FDMH1' 1016 include 'ceedf.h' 1017 1018 parameter (lzq=lx1+1) 1019 parameter (lx=lx1*ly1*lz1,lg=3+3*(ldim-2),lq=lzq**ldim) 1020 common /bpgfactors/ gf(lg*lq,lelt),bmq(lq,lelt),w3mq(lq) 1021 1022 parameter (lt=lx1*ly1*lz1*lelt) 1023 parameter (ld=lxd*lyd*lzd*lelt) 1024 common /vcrns/ u1(lt),r1(lt),r2(lt),r3(lt) 1025 common /vcrny/ e1(lt) 1026 common /vcrvh/ h1(lt),h2(ld),pap(3) 1027 real*8 coords(ldim*lx*lelt) 1028 1029 logical ifh3 1030 integer*8 nnode 1031 integer ceed,err,test 1032 character*64 spec 1033 1034 integer p,q,ncompx,ncompu,enode,lnode 1035 integer vec_p1,vec_ap1,vec_qdata,vec_coords,vec_rhs 1036 integer erstrctu,erstrctx,erstrctw 1037 integer basisu,basisx 1038 integer qf_diffusion,qf_setup 1039 integer op_diffusion,op_setup 1040 integer ii,i,ngeo 1041 real*8 x,y,z 1042 integer*8 offset 1043 1044 external diffusionf,diffsetupf 1045 1046 ifield = 1 1047 nxq = nx1+1 1048 n = nx1*ny1*nz1*nelt 1049 1050 ifsync = .false. 1051 1052C Set up coordinates and mask 1053 ii=0 1054 do j=0,nelt-1 1055 do i=1,lx 1056 ii=ii+1 1057 x = xm1(ii,1,1,1) 1058 y = ym1(ii,1,1,1) 1059 z = zm1(ii,1,1,1) 1060 coords(i+0*lx+3*j*lx)=x 1061 coords(i+1*lx+3*j*lx)=y 1062 coords(i+2*lx+3*j*lx)=z 1063 if ( x.eq.0.or.x.eq.1 1064 $ .or.y.eq.0.or.y.eq.1 1065 $ .or.z.eq.0.or.z.eq.1 ) then 1066 h2(ii)=0. 1067 else 1068 h2(ii)=1. 1069 endif 1070 enddo 1071 enddo 1072 1073C Init ceed library 1074 call get_spec(spec) 1075 call ceedinit(trim(spec)//char(0),ceed,err) 1076 1077 call get_test(test) 1078 1079C Set up Nek geometry data 1080 call geodatq (gf,bmq,w3mq,nxq) ! Set up gf() arrays 1081 1082C Set up true soln 1083 call sin_fld_h1 (e1) 1084 call xmask1 (e1,h2,nelt) 1085 call copy (h1,e1,n) ! Save exact soln in h1 1086 1087C Set up solver parameters 1088 tol = 1e-10 1089 param(22) = tol 1090 maxit = 100 1091 1092 call nekgsync() 1093 1094C Create ceed basis for mesh and computation 1095 p=nx1 1096 q=p+1 1097 ncompu=1 1098 ncompx=ldim 1099 call ceedbasiscreatetensorh1lagrange(ceed,ldim,ncompx,p,q, 1100 $ ceed_gauss,basisx,err) 1101 call ceedbasiscreatetensorh1lagrange(ceed,ldim,ncompu,p,q, 1102 $ ceed_gauss,basisu,err) 1103 1104C Create ceed element restrictions for mesh and computation 1105 enode=p**ldim 1106 lnode=enode*nelt*ncompu 1107 ngeo=(ldim*(ldim+1))/2 1108 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 1109 $ enode,lnode,ldim,erstrctx,err) 1110 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 1111 $ enode,lnode,ncompu,erstrctu,err) 1112 call ceedelemrestrictioncreateidentity(ceed,ceed_notranspose,nelt, 1113 $ q**ldim,nelt*q**ldim,ngeo,erstrctw,err) 1114 1115C Create ceed vectors 1116 call ceedvectorcreate(ceed,lnode,vec_p1,err) 1117 call ceedvectorcreate(ceed,lnode,vec_ap1,err) 1118 call ceedvectorcreate(ceed,lnode,vec_rhs,err) 1119 call ceedvectorcreate(ceed,ldim*lx*nelt,vec_coords,err) 1120 call ceedvectorcreate(ceed,nelt*ngeo*q**ldim,vec_qdata,err) 1121 1122 offset=0 1123 call ceedvectorsetarray(vec_coords,ceed_mem_host, 1124 $ ceed_use_pointer,coords,offset,err) 1125 1126C Create ceed qfunctions for diffsetupf and diffusionf 1127 call ceedqfunctioncreateinterior(ceed,1,diffsetupf, 1128 $ EXAMPLE_DIR 1129 $ //'bps/bps.h:diffsetupf'//char(0),qf_setup,err) 1130 call ceedqfunctionaddinput(qf_setup,'x',ncompx, 1131 $ ceed_eval_interp,err) 1132 call ceedqfunctionaddinput(qf_setup,'dx',ncompx*ldim, 1133 $ ceed_eval_grad,err) 1134 call ceedqfunctionaddinput(qf_setup,'weight',ncompu, 1135 $ ceed_eval_weight,err) 1136 call ceedqfunctionaddoutput(qf_setup,'qdata',ngeo, 1137 $ ceed_eval_none,err) 1138 call ceedqfunctionaddoutput(qf_setup,'rhs',ncompu, 1139 $ ceed_eval_interp,err) 1140 1141 call ceedqfunctioncreateinterior(ceed,1,diffusionf, 1142 $ EXAMPLE_DIR 1143 $ //'bps/bps.h:diffusionf'//char(0),qf_diffusion,err) 1144 call ceedqfunctionaddinput(qf_diffusion,'u',ncompu*ldim, 1145 $ ceed_eval_grad,err) 1146 call ceedqfunctionaddinput(qf_diffusion,'qdata',ngeo, 1147 $ ceed_eval_none,err) 1148 call ceedqfunctionaddoutput(qf_diffusion,'v',ncompu*ldim, 1149 $ ceed_eval_grad,err) 1150 1151C Create ceed operators 1152 call ceedoperatorcreate(ceed,qf_setup, 1153 $ ceed_qfunction_none,ceed_qfunction_none,op_setup,err) 1154 call ceedoperatorsetfield(op_setup,'x',erstrctx, 1155 $ basisx,ceed_vector_active,err) 1156 call ceedoperatorsetfield(op_setup,'dx',erstrctx, 1157 $ basisx,ceed_vector_active,err) 1158 call ceedoperatorsetfield(op_setup,'weight',erstrctx, 1159 $ basisx,ceed_vector_none,err) 1160 call ceedoperatorsetfield(op_setup,'qdata',erstrctw, 1161 $ ceed_basis_collocated,ceed_vector_active,err) 1162 call ceedoperatorsetfield(op_setup,'rhs',erstrctu, 1163 $ basisu,vec_rhs,err) 1164 1165 call ceedoperatorcreate(ceed,qf_diffusion, 1166 $ ceed_qfunction_none,op_diffusion,err) 1167 call ceedoperatorsetfield(op_diffusion,'u',erstrctu, 1168 $ basisu,ceed_vector_active,err) 1169 call ceedoperatorsetfield(op_diffusion,'qdata',erstrctw, 1170 $ ceed_basis_collocated,vec_qdata,err) 1171 call ceedoperatorsetfield(op_diffusion,'v',erstrctu, 1172 $ basisu,ceed_vector_active,err) 1173 1174C Compute setup data 1175 call ceedvectorsetarray(vec_rhs,ceed_mem_host, 1176 $ ceed_use_pointer,r1,offset,err) 1177 call ceedoperatorapply(op_setup,vec_coords,vec_qdata, 1178 $ ceed_request_immediate,err) 1179 1180C Set up true RHS 1181 call dssum (r1,nx1,ny1,nz1) ! r1 1182 call xmask1 (r1,h2,nelt) 1183 1184C Set up algebraic RHS with libCEED 1185 call ceedvectorsetarray(vec_p1,ceed_mem_host, 1186 $ ceed_use_pointer,h1,offset,err) 1187 call ceedvectorsetarray(vec_ap1,ceed_mem_host, 1188 $ ceed_use_pointer,r2,offset,err) 1189 call ceedoperatorapply(op_diffusion,vec_p1,vec_ap1, 1190 $ ceed_request_immediate,err) ! r2 = A_ceed*h1 1191 call dssum (r2,nx1,ny1,nz1) 1192 call xmask1 (r2,h2,nelt) 1193 1194C Set up algebraic RHS with Nek5000 1195 call axhm1 (pap,r3,h1,h1,h2,'bp3') ! r3 = A_nek5k*h1 1196 call dssum (r3,nx1,ny1,nz1) 1197 call xmask1 (r3,h2,nelt) 1198 1199 call nekgsync() 1200 1201C Solve true RHS 1202 if (nid.eq.0) write (6,*) "libCEED true RHS" 1203 tstart = dnekclock() 1204 call cggos(u1,r1,h1,h2,vmult,binvm1,tol,ceed,op_diffusion, 1205 $ vec_p1,vec_ap1,maxit,'bp3') 1206 tstop = dnekclock() 1207 1208C Output 1209 telaps = (tstop-tstart) 1210 maxits = maxit 1211 er1 = glrdif(u1,e1,n) 1212 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 1213 1214 if (test.eq.1.and.nid.eq.0) then 1215 if (maxit>=100) then 1216 write(6,*) "UNCONVERGED CG" 1217 endif 1218 if (dabs(er1)>1e-3) then 1219 write(6,*) "ERROR IS TOO LARGE" 1220 endif 1221 endif 1222 1223 nx = nx1-1 1224 nnode = nelgt ! nnodes 1225 nnode = nnode*(nx**ldim) ! nnodes 1226 nppp = nnode/np ! nnodes/proc 1227 1228 dofps = nnode/telaps ! DOF/sec - scalar form 1229 titers = telaps/maxits ! time per iteration 1230 tppp_s = titers/nppp ! time per iteraton per local point 1231 1232 if (nid.eq.0) write(6,1) 'case scalar:' 1233 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,dofps,titers,tppp_s 1234 1235C Solve libCEED algebraic RHS 1236 if (nid.eq.0) write (6,*) "libCEED algebraic RHS" 1237 maxit = 100 1238 tstart = dnekclock() 1239 call cggos(u1,r2,h1,h2,vmult,binvm1,tol,ceed,op_diffusion, 1240 $ vec_p1,vec_ap1,maxit,'bp3') 1241 tstop = dnekclock() 1242 1243C Output 1244 telaps = (tstop-tstart) 1245 maxits = maxit 1246 er1 = glrdif(u1,e1,n) 1247 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 1248 1249 if (test.eq.1.and.nid.eq.0) then 1250 if (maxit>=100) then 1251 write(6,*) "UNCONVERGED CG" 1252 endif 1253 if (dabs(er1)>1e-9) then 1254 write(6,*) "ERROR IS TOO LARGE" 1255 endif 1256 endif 1257 1258 nx = nx1-1 1259 nnode = nelgt ! nnodes 1260 nnode = nnode*(nx**ldim) ! nnodes 1261 nppp = nnode/np ! nnodes/proc 1262 1263 dofps = nnode/telaps ! DOF/sec - scalar form 1264 titers = telaps/maxits ! time per iteration 1265 tppp_s = titers/nppp ! time per iteraton per local point 1266 1267 if (nid.eq.0) write(6,1) 'case scalar:' 1268 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,dofps,titers,tppp_s 1269 1270C Solve Nek5000 algebraic RHS 1271 if (nid.eq.0) write (6,*) "Nek5000 algebraic RHS" 1272 maxit = 100 1273 tstart = dnekclock() 1274 call cggos(u1,r3,h1,h2,vmult,binvm1,tol,ceed,op_diffusion, 1275 $ vec_p1,vec_ap1,maxit,'bp3') 1276 tstop = dnekclock() 1277 1278C Output 1279 telaps = (tstop-tstart) 1280 maxits = maxit 1281 er1 = glrdif(u1,e1,n) 1282 if (nid.eq.0) write(6,3) lx1,nelgv,er1,' error ',maxit 1283 1284 if (test.eq.1.and.nid.eq.0) then 1285 if (maxit>=100) then 1286 write(6,*) "UNCONVERGED CG" 1287 endif 1288 if (dabs(er1)>1e-9) then 1289 write(6,*) "ERROR IS TOO LARGE" 1290 endif 1291 endif 1292 1293 nx = nx1-1 1294 nnode = nelgt ! nnodes 1295 nnode = nnode*(nx**ldim) ! nnodes 1296 nppp = nnode/np ! nnodes/proc 1297 1298 nodepss = nnode/telaps ! DOF/sec - scalar form 1299 titers = telaps/maxits ! time per iteration 1300 tppp_s = titers/nppp ! time per iteraton per local point 1301 1302 if (nid.eq.0) write(6,1) 'case scalar:' 1303 $ ,np,nx,nelt,nelgt,nnode,nppp,maxits,telaps,nodepss,titers,tppp_s 1304 1305 1 format(a12,i7,i3,i7,i10,i14,i10,i4,1p4e13.5) 1306 3 format(i3,i9,e12.4,1x,a8,i9) 1307 1308C Destroy ceed handles 1309 call ceedvectordestroy(vec_p1,err) 1310 call ceedvectordestroy(vec_ap1,err) 1311 call ceedvectordestroy(vec_rhs,err) 1312 call ceedvectordestroy(vec_qdata,err) 1313 call ceedvectordestroy(vec_coords,err) 1314 call ceedelemrestrictiondestroy(erstrctu,err) 1315 call ceedelemrestrictiondestroy(erstrctx,err) 1316 call ceedelemrestrictiondestroy(erstrctw,err) 1317 call ceedbasisdestroy(basisu,err) 1318 call ceedbasisdestroy(basisx,err) 1319 call ceedqfunctiondestroy(qf_setup,err) 1320 call ceedqfunctiondestroy(qf_diffusion,err) 1321 call ceedoperatordestroy(op_setup,err) 1322 call ceedoperatordestroy(op_diffusion,err) 1323 call ceeddestroy(ceed,err) 1324 1325 return 1326 end 1327C----------------------------------------------------------------------- 1328 subroutine cggos(u1,r1,h1,h2,rmult,binv,tin,ceed,ceed_op,vec_p1, 1329 $ vec_ap1,maxit,bpname) 1330C Scalar conjugate gradient iteration for solution of uncoupled 1331C Helmholtz equations 1332C Input: r1,h1,h2,rmult,binv,tin,ceed,ceed_op,vec_p1,vec_ap1,bpname 1333C Output: u1,maxit 1334 include 'SIZE' 1335 include 'TOTAL' 1336 include 'DOMAIN' 1337 include 'FDMH1' 1338 character*3 bpname 1339 1340C INPUT: rhs1 - rhs 1341C h1 - exact solution 1342 1343 parameter (lt=lx1*ly1*lz1*lelt) 1344 parameter (ld=lxd*lyd*lzd*lelt) 1345 real*8 u1(lt),r1(lt),h1(lt),h2(lt) 1346 real*8 rmult(1),binv(1) 1347 integer ceed,ceed_op,vec_ap1,vec_p1 1348 common /scrcg/ dpc(lt),p1(lt),z1(lt) 1349 common /scrca/ wv(4),wk(4),rpp1(4),rpp2(4),alph(4),beta(4),pap(4) 1350 1351 real*8 ap1(lt) 1352 equivalence (ap1,z1) 1353 1354 vol = volfld(ifield) 1355 nel = nelfld(ifield) 1356 nxyz = lx1*ly1*lz1 1357 n = nxyz*nel 1358 nx = nx1-1 ! Polynomial order (just for i/o) 1359 1360 tol=tin 1361 1362 if(bpname.ne.'bp1') then 1363 call setprecn_bp3(dpc,h1,h2) ! Set up diagional pre-conidtioner 1364 else 1365 call setprecn_bp1(dpc,h1,h2) ! Set up diagional pre-conidtioner 1366 endif 1367 1368 call rzero (u1,n) ! Initialize solution 1369 1370 wv(1)=0 1371 do i=1,n 1372 s=rmult(i) ! -1 1373 p1(i)=dpc(i)*r1(i) ! p = M r T 1374 wv(1)=wv(1)+s*p1(i)*r1(i) ! r p 1375 enddo 1376 call gop(wv(1),wk,'+ ',1) 1377 rpp1(1) = wv (1) 1378 1379 do 1000 iter=1,maxit 1380 call axhm1_ceed (pap,ap1,p1,h1,h2,ceed,ceed_op, 1381 $ vec_ap1,vec_p1) 1382 call dssum (ap1,nx1,ny1,nz1) 1383 if (bpname.ne.'bp1') call xmask1(ap1,h2,nel) 1384 1385 call gop (pap,wk,'+ ',1) 1386 alph(1) = rpp1(1)/pap(1) 1387 1388 do i=1,n 1389 u1(i)=u1(i)+alph(1)* p1(i) 1390 r1(i)=r1(i)-alph(1)*ap1(i) 1391 enddo 1392 1393C tolerance check here 1394 call rzero(wv,2) 1395 do i=1,n 1396 wv(1)=wv(1)+r1(i)*r1(i) ! L2 error estimate 1397 z1(i)=dpc(i)*r1(i) ! z = M r 1398 wv(2)=wv(2)+rmult(i)*z1(i)*r1(i) ! r z 1399 enddo 1400 call gop(wv,wk,'+ ',2) 1401 1402C if (nio.eq.0) write(6,1) ifield,istep,iter,nx,(wv(k),k=1,1) 1403 1 format(i2,i9,i5,i4,1p1e12.4,' cggos') 1404 1405 enorm=sqrt(wv(1)) 1406 if (enorm.lt.tol) then 1407 ifin = iter 1408 if (nio.eq.0) write(6,3000) istep,ifin,enorm,tol 1409 goto 9999 1410 endif 1411C if (nio.eq.0) write(6,2) iter,enorm,alph(1),pap(1),'alpha' 1412 2 format(i5,1p3e12.4,2x,a5) 1413 1414 rpp2(1)=rpp1(1) 1415 rpp1(1)=wv (2) 1416 beta1 =rpp1(1)/rpp2(1) 1417 do i=1,n 1418 p1(i)=z1(i) + beta1*p1(i) 1419 enddo 1420 1421 1000 continue 1422 1423 rbnorm=sqrt(wv(1)) 1424 if (nio.eq.0) write (6,3001) istep,iter,rbnorm,tol 1425 iter = iter-1 1426 1427 9999 continue 1428 1429 maxit=iter 1430 1431 3000 format(i12,1x,'cggo scalar:',i6,1p5e13.4) 1432 3001 format(2i6,' Unconverged cggo scalar: rbnorm =',1p2e13.6) 1433 1434 return 1435 end 1436C----------------------------------------------------------------------- 1437 subroutine axhm1_ceed(pap,ap1,p1,h1,h2,ceed,ceed_op, 1438 $ vec_ap1,vec_p1) 1439C Vector conjugate gradient matvec for solution of uncoupled 1440C Helmholtz equations 1441C Input: pap,p1,h1,h2,bpname,ceed,ceed_op,vec_ap1,vec_p1 1442C Output: ap1 1443 include 'SIZE' 1444 include 'TOTAL' 1445 include 'ceedf.h' 1446 1447 parameter (lx=lx1*ly1*lz1,lg=3+3*(ldim-2)) 1448 real*8 gf(lg,lx,lelt) ! Equivalence new gf() data 1449 equivalence (gf,g1m1) ! layout to g1m1...g6m1 1450 1451 real*8 pap(3) 1452 real*8 ap1(lx,lelt) 1453 real*8 p1(lx,lelt) 1454 real*8 h1(lx,lelt),h2(lx,lelt) 1455 integer ceed,ceed_op,vec_ap1,vec_p1,err 1456 integer i,e 1457 integer*8 offset 1458 1459 offset=0 1460 call ceedvectorsetarray(vec_p1,ceed_mem_host,ceed_use_pointer, 1461 $ p1,offset,err) 1462 call ceedvectorsetarray(vec_ap1,ceed_mem_host,ceed_use_pointer, 1463 $ ap1,offset,err) 1464 1465 call ceedoperatorapply(ceed_op,vec_p1,vec_ap1, 1466 $ ceed_request_immediate,err) 1467 1468 call ceedvectorsyncarray(vec_ap1,ceed_mem_host,err) 1469 1470 pap(1)=0. 1471 1472 do e=1,nelt 1473 do i=1,lx 1474 pap(1)=pap(1)+ap1(i,e)*p1(i,e) 1475 enddo 1476 enddo 1477 1478 return 1479 end 1480C----------------------------------------------------------------------- 1481 subroutine ax_e_bp1(w,u,g,h1,h2,b,ju,us,ut) 1482C Local matrix-vector for solution of BP3 (stiffness matrix) 1483C Input: u,g,h1,h2,b,ju,us,ut Output: w 1484 include 'SIZE' 1485 include 'TOTAL' 1486 1487 parameter (lxyz=lx1*ly1*lz1,lg=3+3*(ldim-2)) 1488 real*8 w(lxyz),u(lxyz),g(lg,lxyz),h1(lxyz),h2(lxyz),b(lxyz) 1489 real*8 ju(lxyz),us(lxyz),ut(lxyz) 1490 1491 nxq = nx1+1 ! Number of quadrature points 1492 1493 lxyzq = nxq**ldim 1494 1495 call intp_rstd (ju,u,lx1,nxq,if3d,0) ! 0 --> Fwd interpolation 1496 do i=1,lxyzq 1497 ju(i)=ju(i)*h2(i) !! h2 must be on the fine grid, w/ quad wts 1498 enddo 1499 call intp_rstd (w,ju,lx1,nxq,if3d,1) ! 1 --> ju-->u 1500 1501 return 1502 end 1503C----------------------------------------------------------------------- 1504 subroutine axhm1_bp1(pap,ap1,p1,h1,h2) 1505C Vector conjugate gradient matvec for solution of BP1 (mass matrix) 1506C Input: pap,p1,h1,h2 Output: ap1 1507 include 'SIZE' 1508 include 'TOTAL' 1509 1510 parameter (lx=lx1*ly1*lz1,lg=3+3*(ldim-2)) 1511 real*8 gf(lg,lx,lelt) ! Equivalence new gf() data 1512 equivalence (gf,g1m1) ! layout to g1m1...g6m1 1513 1514 real*8 pap(3) 1515 real*8 ap1(lx,lelt) 1516 real*8 p1(lx,lelt) 1517 real*8 h1(lx,lelt),h2(lx,lelt) 1518 1519 real*8 ur(lx),us(lx),ut(lx) 1520 common /ctmp1/ ur,us,ut 1521 1522 integer e 1523 1524 pap(1)=0. 1525 1526 k=1 1527 nxq = nx1+1 1528 1529 do e=1,nelt 1530 1531 call ax_e_bp1(ap1(1,e),p1(1,e),gf(1,1,e),h1(1,e),h2(k,1) 1532 $ ,bm1(1,1,1,e),ur,us,ut) 1533 do i=1,lx 1534 pap(1)=pap(1)+ap1(i,e)*p1(i,e) 1535 enddo 1536 k=k+nxq*nxq*nxq 1537 1538 enddo 1539 1540 return 1541 end 1542C----------------------------------------------------------------------- 1543 subroutine ax_e_bp3(w,u,g,ur,us,ut,wk) 1544C Local matrix-vector for solution of BP3 (stiffness matrix) 1545C Input: u,g,ur,us,ut,wk Output: w 1546 include 'SIZE' 1547 include 'TOTAL' 1548 1549 parameter (lzq=lx1+1,lxyz=lx1*lx1*lx1,lxyzq=lzq*lzq*lzq) 1550 1551 common /ctmp0/ tmp(lxyzq) 1552 common /dxmfine/ dxmq(lzq,lzq),dxtmq(lzq,lzq) 1553 1554 real*8 ur(lxyzq),us(lxyzq),ut(lxyzq),wk(lxyzq) 1555 real*8 w(lxyz),u(lxyz),g(2*ldim,lxyzq) 1556 1557 n = lzq-1 1558 1559 call intp_rstd (wk,u,lx1,lzq,if3d,0) ! 0 --> Fwd interpolation 1560 call loc_grad3 (ur,us,ut,wk,n,dxmq,dxtmq) 1561 1562 do i=1,lxyzq 1563 wr = g(1,i)*ur(i) + g(2,i)*us(i) + g(3,i)*ut(i) 1564 ws = g(2,i)*ur(i) + g(4,i)*us(i) + g(5,i)*ut(i) 1565 wt = g(3,i)*ur(i) + g(5,i)*us(i) + g(6,i)*ut(i) 1566 ur(i) = wr 1567 us(i) = ws 1568 ut(i) = wt 1569 enddo 1570 1571 call loc_grad3t (wk,ur,us,ut,n,dxmq,dxtmq,tmp) 1572 call intp_rstd (w,wk,lx1,lzq,if3d,1) ! 1 --> ju-->u 1573 1574 return 1575 end 1576C----------------------------------------------------------------------- 1577 subroutine axhm1_bp3(pap,ap1,p1,h1,h2) 1578C Vector conjugate gradient matvec for solution of BP3 (stiffness matrix) 1579C Input: pap,p1,h1,h2 Output: ap1 1580 include 'SIZE' 1581 include 'TOTAL' 1582 1583 parameter (lzq=lx1+1) 1584 parameter (lx=lx1*ly1*lz1,lg=3+3*(ldim-2),lq=lzq**ldim) 1585 common /bpgfactors/ gf(lg,lq,lelt),bmq(lq,lelt),w3mq(lq) 1586 1587 real*8 pap(3) 1588 real*8 ap1(lx,lelt) 1589 real*8 p1(lx,lelt) 1590 real*8 h1(lx,lelt),h2(lx,lelt) 1591 1592 common /ctmp1/ ur,us,ut,wk 1593 real*8 ur(lq),us(lq),ut(lq),wk(lq) 1594 1595 integer e 1596 1597 pap(1)=0. 1598 1599 do e=1,nelt 1600 1601 call ax_e_bp3(ap1(1,e),p1(1,e),gf(1,1,e),ur,us,ut,wk) 1602 do i=1,lx 1603 pap(1)=pap(1)+p1(i,e)*ap1(i,e) 1604 enddo 1605 1606 enddo 1607 1608 return 1609 end 1610C----------------------------------------------------------------------- 1611 subroutine axhm1(pap,ap1,p1,h1,h2,bpname) 1612C Vector conjugate gradient matvec for solution of uncoupled 1613C Helmholtz equations 1614C Input: pap,p1,h1,h2,bpname Output: ap1 1615 include 'SIZE' 1616 include 'TOTAL' 1617 1618 parameter (lx=lx1*ly1*lz1) 1619 1620 real*8 pap(3),ap1(lx,lelt),p1(lx,lelt) 1621 real*8 h1(lx,lelt),h2(lx,lelt) 1622 1623 character*3 bpname 1624 1625 if (bpname.eq.'bp1') then 1626 call axhm1_bp1(pap,ap1,p1,h1,h2) 1627 1628 elseif (bpname.eq.'bp3') then 1629 call axhm1_bp3(pap,ap1,p1,h1,h2) 1630 1631 else 1632 write(6,*) bpname,' axhm1 bpname error' 1633 stop 1634 1635 endif 1636 1637 return 1638 end 1639C----------------------------------------------------------------------- 1640 subroutine get_bp(bp) 1641C Get BP to run 1642C Input: Output: bp 1643 integer i,bp 1644 character*64 bpval 1645 1646 bp=0 1647 if(iargc().ge.1) then 1648 call getarg(1,bpval) 1649 endif 1650 if(bpval.eq."bp1") then 1651 bp=1 1652 elseif(bpval.eq."bp3") then 1653 bp=3 1654 endif 1655 1656 return 1657 end 1658C----------------------------------------------------------------------- 1659 subroutine get_spec(spec) 1660C Get CEED backend specification 1661C Input: Output: spec 1662 integer i 1663 character*64 spec 1664 1665 spec = '/cpu/self' 1666 if(iargc().ge.2) then 1667 call getarg(2,spec) 1668 endif 1669 1670 return 1671 end 1672C----------------------------------------------------------------------- 1673 subroutine get_test(test) 1674C Get test mode flag 1675C Input: Output: test 1676 integer i,test 1677 character*64 testval 1678 1679 test=0 1680 if(iargc().ge.3) then 1681 call getarg(3,testval) 1682 endif 1683 if(testval.eq."test") then 1684 test=1 1685 endif 1686 1687 return 1688 end 1689C----------------------------------------------------------------------- 1690