| 1 | c-----------------------------------------------------------------------
|
|---|
| 2 | subroutine nek_init(comm)
|
|---|
| 3 | c
|
|---|
| 4 | include 'SIZE'
|
|---|
| 5 | include 'TOTAL'
|
|---|
| 6 | include 'DOMAIN'
|
|---|
| 7 | c
|
|---|
| 8 | include 'OPCTR'
|
|---|
| 9 | include 'CTIMER'
|
|---|
| 10 |
|
|---|
| 11 | C used scratch arrays
|
|---|
| 12 | C NOTE: no initial declaration needed. Linker will take
|
|---|
| 13 | c care about the size of the CBs automatically
|
|---|
| 14 | c
|
|---|
| 15 | c COMMON /CTMP1/ DUMMY1(LCTMP1)
|
|---|
| 16 | c COMMON /CTMP0/ DUMMY0(LCTMP0)
|
|---|
| 17 | c
|
|---|
| 18 | c COMMON /SCRNS/ DUMMY2(LX1,LY1,LZ1,LELT,7)
|
|---|
| 19 | c COMMON /SCRUZ/ DUMMY3(LX1,LY1,LZ1,LELT,4)
|
|---|
| 20 | c COMMON /SCREV/ DUMMY4(LX1,LY1,LZ1,LELT,2)
|
|---|
| 21 | c COMMON /SCRVH/ DUMMY5(LX1,LY1,LZ1,LELT,2)
|
|---|
| 22 | c COMMON /SCRMG/ DUMMY6(LX1,LY1,LZ1,LELT,4)
|
|---|
| 23 | c COMMON /SCRCH/ DUMMY7(LX1,LY1,LZ1,LELT,2)
|
|---|
| 24 | c COMMON /SCRSF/ DUMMY8(LX1,LY1,LZ1,LELT,3)
|
|---|
| 25 | c COMMON /SCRCG/ DUMM10(LX1,LY1,LZ1,LELT,1)
|
|---|
| 26 |
|
|---|
| 27 | integer comm
|
|---|
| 28 | common /nekmpi/ mid,mp,nekcomm,nekgroup,nekreal
|
|---|
| 29 |
|
|---|
| 30 | common /rdump/ ntdump
|
|---|
| 31 |
|
|---|
| 32 | real kwave2
|
|---|
| 33 | logical ifemati
|
|---|
| 34 |
|
|---|
| 35 | real rtest
|
|---|
| 36 | integer itest
|
|---|
| 37 | integer*8 itest8
|
|---|
| 38 | character ctest
|
|---|
| 39 | logical ltest
|
|---|
| 40 |
|
|---|
| 41 | common /c_is1/ glo_num(lx1 * ly1 * lz1, lelt)
|
|---|
| 42 | common /ivrtx/ vertex((2 ** ldim) * lelt)
|
|---|
| 43 | integer*8 glo_num, ngv
|
|---|
| 44 | integer vertex
|
|---|
| 45 |
|
|---|
| 46 | ! set word size for REAL
|
|---|
| 47 | wdsize = sizeof(rtest)
|
|---|
| 48 | ! set word size for INTEGER
|
|---|
| 49 | isize = sizeof(itest)
|
|---|
| 50 | ! set word size for INTEGER*8
|
|---|
| 51 | isize8 = sizeof(itest8)
|
|---|
| 52 | ! set word size for LOGICAL
|
|---|
| 53 | lsize = sizeof(ltest)
|
|---|
| 54 | ! set word size for CHARACTER
|
|---|
| 55 | csize = sizeof(ctest)
|
|---|
| 56 |
|
|---|
| 57 | call setupcomm(comm,newcomm,newcommg,'','')
|
|---|
| 58 | intracomm = newcomm ! within a session
|
|---|
| 59 | nekcomm = newcomm
|
|---|
| 60 | iglobalcomm = newcommg ! across all sessions
|
|---|
| 61 | call iniproc()
|
|---|
| 62 |
|
|---|
| 63 | if (nid.eq.nio) call printHeader
|
|---|
| 64 |
|
|---|
| 65 | etimes = dnekclock()
|
|---|
| 66 | istep = 0
|
|---|
| 67 |
|
|---|
| 68 | call opcount(1)
|
|---|
| 69 |
|
|---|
| 70 | call initdim ! Initialize / set default values.
|
|---|
| 71 | call initdat
|
|---|
| 72 | call files
|
|---|
| 73 |
|
|---|
| 74 | call readat ! Read .rea +map file
|
|---|
| 75 |
|
|---|
| 76 | call setvar ! Initialize most variables
|
|---|
| 77 |
|
|---|
| 78 | instep=1 ! Check for zero steps
|
|---|
| 79 | if (nsteps.eq.0 .and. fintim.eq.0.) instep=0
|
|---|
| 80 |
|
|---|
| 81 | igeom = 2
|
|---|
| 82 | call setup_topo ! Setup domain topology
|
|---|
| 83 |
|
|---|
| 84 | call genwz ! Compute GLL points, weights, etc.
|
|---|
| 85 |
|
|---|
| 86 | if(nio.eq.0) write(6,*) 'call usrdat'
|
|---|
| 87 | call usrdat
|
|---|
| 88 | if(nio.eq.0) write(6,'(A,/)') ' done :: usrdat'
|
|---|
| 89 |
|
|---|
| 90 | call gengeom(igeom) ! Generate geometry, after usrdat
|
|---|
| 91 |
|
|---|
| 92 | if (ifmvbd) call setup_mesh_dssum ! Set mesh dssum (needs geom)
|
|---|
| 93 |
|
|---|
| 94 | if(nio.eq.0) write(6,*) 'call usrdat2'
|
|---|
| 95 | call usrdat2
|
|---|
| 96 | if(nio.eq.0) write(6,'(A,/)') ' done :: usrdat2'
|
|---|
| 97 |
|
|---|
| 98 | call fix_geom
|
|---|
| 99 | call geom_reset(1) ! recompute Jacobians, etc.
|
|---|
| 100 |
|
|---|
| 101 | call vrdsmsh ! verify mesh topology
|
|---|
| 102 | call mesh_metrics ! print some metrics
|
|---|
| 103 |
|
|---|
| 104 | call setlog(.true.) ! Initalize logical flags
|
|---|
| 105 |
|
|---|
| 106 | if (ifneknekc) call neknek_setup
|
|---|
| 107 |
|
|---|
| 108 | call bcmask ! Set BC masks for Dirichlet boundaries.
|
|---|
| 109 |
|
|---|
| 110 | if (fintim.ne.0.0 .or. nsteps.ne.0)
|
|---|
| 111 | $ call geneig(igeom) ! eigvals for tolerances
|
|---|
| 112 |
|
|---|
| 113 | call dg_setup ! Setup DG, if dg flag is set.
|
|---|
| 114 |
|
|---|
| 115 | if (ifflow.and.iftran) then ! Init pressure solver
|
|---|
| 116 | if (fintim.ne.0 .or. nsteps.ne.0) call prinit
|
|---|
| 117 | endif
|
|---|
| 118 |
|
|---|
| 119 | if(ifcvode) call cv_setsize
|
|---|
| 120 |
|
|---|
| 121 | if(nio.eq.0) write(6,*) 'call usrdat3'
|
|---|
| 122 | call usrdat3
|
|---|
| 123 | if(nio.eq.0) write(6,'(A,/)') ' done :: usrdat3'
|
|---|
| 124 |
|
|---|
| 125 | #ifdef CMTNEK
|
|---|
| 126 | call nek_cmt_init
|
|---|
| 127 | #endif
|
|---|
| 128 |
|
|---|
| 129 | call setics
|
|---|
| 130 | call setprop
|
|---|
| 131 |
|
|---|
| 132 | if (instep.ne.0) then
|
|---|
| 133 | if (ifneknekc) call neknek_exchange
|
|---|
| 134 | if (ifneknekc) call chk_outflow
|
|---|
| 135 |
|
|---|
| 136 | if (nio.eq.0) write(6,*) 'call userchk'
|
|---|
| 137 | call userchk
|
|---|
| 138 | if(nio.eq.0) write(6,'(A,/)') ' done :: userchk'
|
|---|
| 139 | endif
|
|---|
| 140 |
|
|---|
| 141 | call setprop ! call again because input has changed in userchk
|
|---|
| 142 |
|
|---|
| 143 | if (ifcvode .and. nsteps.gt.0) call cv_init
|
|---|
| 144 |
|
|---|
| 145 | call comment
|
|---|
| 146 | call sstest (isss)
|
|---|
| 147 |
|
|---|
| 148 | call dofcnt
|
|---|
| 149 |
|
|---|
| 150 | jp = 0 ! Set perturbation field count to 0 for baseline flow
|
|---|
| 151 | p0thn = p0th
|
|---|
| 152 |
|
|---|
| 153 | call in_situ_init()
|
|---|
| 154 |
|
|---|
| 155 | call time00 ! Initalize timers to ZERO
|
|---|
| 156 | call opcount(2)
|
|---|
| 157 |
|
|---|
| 158 | ntdump=0
|
|---|
| 159 | if (timeio.ne.0.0) ntdump = int( time/timeio )
|
|---|
| 160 |
|
|---|
| 161 | tinit = dnekclock_sync() - etimes
|
|---|
| 162 | if (nio.eq.0) then
|
|---|
| 163 | write (6,*) ' '
|
|---|
| 164 | if (time.ne.0.0) write (6,'(a,e14.7)') ' Initial time:',time
|
|---|
| 165 | write (6,'(a,g13.5,a)')
|
|---|
| 166 | & ' Initialization successfully completed ', tinit, ' sec'
|
|---|
| 167 | endif
|
|---|
| 168 |
|
|---|
| 169 | return
|
|---|
| 170 | end
|
|---|
| 171 | c-----------------------------------------------------------------------
|
|---|
| 172 | subroutine nek_solve
|
|---|
| 173 |
|
|---|
| 174 | include 'SIZE'
|
|---|
| 175 | include 'TSTEP'
|
|---|
| 176 | include 'INPUT'
|
|---|
| 177 | include 'CTIMER'
|
|---|
| 178 |
|
|---|
| 179 | call nekgsync()
|
|---|
| 180 |
|
|---|
| 181 | if (instep.eq.0) then
|
|---|
| 182 | if(nid.eq.0) write(6,'(/,A,/,A,/)')
|
|---|
| 183 | & ' nsteps=0 -> skip time loop',
|
|---|
| 184 | & ' running solver in post processing mode'
|
|---|
| 185 | else
|
|---|
| 186 | if(nio.eq.0) write(6,'(/,A,/)') 'Starting time loop ...'
|
|---|
| 187 | endif
|
|---|
| 188 |
|
|---|
| 189 | isyc = 0
|
|---|
| 190 | if(ifsync) isyc=1
|
|---|
| 191 | itime = 0
|
|---|
| 192 | #ifdef TIMER
|
|---|
| 193 | itime = 1
|
|---|
| 194 | #endif
|
|---|
| 195 |
|
|---|
| 196 | ! start measurements
|
|---|
| 197 | dtmp = dnekgflops()
|
|---|
| 198 |
|
|---|
| 199 | istep = 0
|
|---|
| 200 | msteps = 1
|
|---|
| 201 |
|
|---|
| 202 | irstat = int(param(120))
|
|---|
| 203 |
|
|---|
| 204 | do kstep=1,nsteps,msteps
|
|---|
| 205 | call nek__multi_advance(kstep,msteps)
|
|---|
| 206 | if(kstep.ge.nsteps) lastep = 1
|
|---|
| 207 | call check_ioinfo
|
|---|
| 208 | call set_outfld
|
|---|
| 209 | etime1 = dnekclock()
|
|---|
| 210 | call userchk
|
|---|
| 211 | tuchk = tuchk + dnekclock()-etime1
|
|---|
| 212 | call prepost (ifoutfld,'his')
|
|---|
| 213 | call in_situ_check()
|
|---|
| 214 | if (mod(kstep,irstat).eq.0 .and. lastep.eq.0) call runstat
|
|---|
| 215 | if (lastep .eq. 1) goto 1001
|
|---|
| 216 | enddo
|
|---|
| 217 | 1001 lastep=1
|
|---|
| 218 |
|
|---|
| 219 | call comment
|
|---|
| 220 |
|
|---|
| 221 | c check for post-processing mode
|
|---|
| 222 | if (instep.eq.0) then
|
|---|
| 223 | nsteps=0
|
|---|
| 224 | istep=0
|
|---|
| 225 | if(nio.eq.0) write(6,*) 'call userchk'
|
|---|
| 226 | call userchk
|
|---|
| 227 | if(nio.eq.0) write(6,*) 'done :: userchk'
|
|---|
| 228 | call prepost (.true.,'his')
|
|---|
| 229 | else
|
|---|
| 230 | if (nio.eq.0) write(6,'(/,A,/)')
|
|---|
| 231 | $ 'end of time-step loop'
|
|---|
| 232 | endif
|
|---|
| 233 |
|
|---|
| 234 |
|
|---|
| 235 | RETURN
|
|---|
| 236 | END
|
|---|
| 237 |
|
|---|
| 238 | c-----------------------------------------------------------------------
|
|---|
| 239 | subroutine nek_advance
|
|---|
| 240 |
|
|---|
| 241 | include 'SIZE'
|
|---|
| 242 | include 'TOTAL'
|
|---|
| 243 | include 'CTIMER'
|
|---|
| 244 |
|
|---|
| 245 | common /cgeom/ igeom
|
|---|
| 246 |
|
|---|
| 247 | ntot = lx1*ly1*lz1*nelv
|
|---|
| 248 |
|
|---|
| 249 | call nekgsync
|
|---|
| 250 |
|
|---|
| 251 | call setup_convect(2) ! Save conv vel
|
|---|
| 252 |
|
|---|
| 253 | if (iftran) call settime
|
|---|
| 254 | if (ifmhd ) call cfl_check
|
|---|
| 255 | call setsolv
|
|---|
| 256 | call comment
|
|---|
| 257 |
|
|---|
| 258 | #ifdef CMTNEK
|
|---|
| 259 | if (nio.eq.0.and.istep.le.1) write(6,*) 'CMT branch active'
|
|---|
| 260 | call cmt_nek_advance
|
|---|
| 261 | return
|
|---|
| 262 | #endif
|
|---|
| 263 |
|
|---|
| 264 | if (ifsplit) then ! PN/PN formulation
|
|---|
| 265 |
|
|---|
| 266 | do igeom=1,ngeom
|
|---|
| 267 |
|
|---|
| 268 | if (ifneknekc .and. igeom.gt.2) then
|
|---|
| 269 | if (ifneknekm.and.igeom.eq.3) call neknek_setup
|
|---|
| 270 | call neknek_exchange
|
|---|
| 271 | endif
|
|---|
| 272 |
|
|---|
| 273 | ! call here before we overwrite wx
|
|---|
| 274 | if (ifheat .and. ifcvode) call heat_cvode (igeom)
|
|---|
| 275 |
|
|---|
| 276 | if (ifgeom) then
|
|---|
| 277 | call gengeom (igeom)
|
|---|
| 278 | call geneig (igeom)
|
|---|
| 279 | endif
|
|---|
| 280 |
|
|---|
| 281 | if (ifheat) call heat (igeom)
|
|---|
| 282 |
|
|---|
| 283 | if (igeom.eq.2) then
|
|---|
| 284 | call setprop
|
|---|
| 285 | call rzero(qtl,ntot)
|
|---|
| 286 | if (iflomach) call qthermal
|
|---|
| 287 | endif
|
|---|
| 288 |
|
|---|
| 289 | if (ifflow) call fluid (igeom)
|
|---|
| 290 | if (ifmvbd) call meshv (igeom)
|
|---|
| 291 | if (igeom.eq.ngeom.and.filterType.eq.1)
|
|---|
| 292 | $ call q_filter(param(103))
|
|---|
| 293 |
|
|---|
| 294 | enddo
|
|---|
| 295 |
|
|---|
| 296 | else ! PN-2/PN-2 formulation
|
|---|
| 297 | call setprop
|
|---|
| 298 | do igeom=1,ngeom
|
|---|
| 299 |
|
|---|
| 300 | if (ifneknekc .and. igeom.gt.2) then
|
|---|
| 301 | if (ifneknekm.and.igeom.eq.3) call neknek_setup
|
|---|
| 302 | call neknek_exchange
|
|---|
| 303 | endif
|
|---|
| 304 |
|
|---|
| 305 | ! call here before we overwrite wx
|
|---|
| 306 | if (ifheat .and. ifcvode) call heat_cvode (igeom)
|
|---|
| 307 |
|
|---|
| 308 | if (ifgeom) then
|
|---|
| 309 | if (.not.ifrich) call gengeom (igeom)
|
|---|
| 310 | call geneig (igeom)
|
|---|
| 311 | endif
|
|---|
| 312 |
|
|---|
| 313 | if (ifmhd) then
|
|---|
| 314 | if (ifheat) call heat (igeom)
|
|---|
| 315 | call induct (igeom)
|
|---|
| 316 | elseif (ifpert) then
|
|---|
| 317 | if (ifbase.and.ifheat) call heat (igeom)
|
|---|
| 318 | if (ifbase.and.ifflow) call fluid (igeom)
|
|---|
| 319 | if (ifflow) call fluidp (igeom)
|
|---|
| 320 | if (ifheat) call heatp (igeom)
|
|---|
| 321 | else ! std. nek case
|
|---|
| 322 | if (ifheat) call heat (igeom)
|
|---|
| 323 | if (ifflow) call fluid (igeom)
|
|---|
| 324 | if (ifmvbd) call meshv (igeom)
|
|---|
| 325 | endif
|
|---|
| 326 | if (igeom.eq.ngeom.and.filterType.eq.1)
|
|---|
| 327 | $ call q_filter(param(103))
|
|---|
| 328 | enddo
|
|---|
| 329 | endif
|
|---|
| 330 |
|
|---|
| 331 | return
|
|---|
| 332 | end
|
|---|
| 333 |
|
|---|
| 334 | c-----------------------------------------------------------------------
|
|---|
| 335 | subroutine nek_end
|
|---|
| 336 |
|
|---|
| 337 | include 'SIZE'
|
|---|
| 338 | include 'TOTAL'
|
|---|
| 339 |
|
|---|
| 340 | if(instep.ne.0) call runstat
|
|---|
| 341 |
|
|---|
| 342 | c if (ifstrs) then
|
|---|
| 343 | c call fgslib_crs_free(xxth_strs)
|
|---|
| 344 | c else
|
|---|
| 345 | c call fgslib_crs_free(xxth(1))
|
|---|
| 346 | c endif
|
|---|
| 347 |
|
|---|
| 348 | call in_situ_end()
|
|---|
| 349 | call exitt0()
|
|---|
| 350 |
|
|---|
| 351 | return
|
|---|
| 352 | end
|
|---|
| 353 | c-----------------------------------------------------------------------
|
|---|
| 354 | subroutine nek__multi_advance(kstep,msteps)
|
|---|
| 355 |
|
|---|
| 356 | include 'SIZE'
|
|---|
| 357 | include 'TOTAL'
|
|---|
| 358 |
|
|---|
| 359 | do i=1,msteps
|
|---|
| 360 | istep = istep+i
|
|---|
| 361 | call nek_advance
|
|---|
| 362 |
|
|---|
| 363 | if (ifneknekc) then
|
|---|
| 364 | call neknek_exchange
|
|---|
| 365 | call bcopy
|
|---|
| 366 | call chk_outflow
|
|---|
| 367 | endif
|
|---|
| 368 | enddo
|
|---|
| 369 |
|
|---|
| 370 | return
|
|---|
| 371 | end
|
|---|
| 372 | c-----------------------------------------------------------------------
|
|---|