5 subroutine whalo1(level, start, end, commPressure, commGamma, &
24 integer(kind=intType),
intent(in) :: level, start, end
25 logical,
intent(in) :: commPressure, commGamma, commViscous
29 integer(kind=intType) :: nn, mm, ll
31 logical :: correctForK, commLamVis, commEddyVis, commVarGamma
36 if (
viscous .and. commviscous) commlamvis = .true.
39 if (
eddymodel .and. commviscous) commeddyvis = .true.
43 commvargamma = .false.
49 call whalo1to1(level, start,
end, commPressure, commVarGamma, &
50 commlamvis, commeddyvis, commpatterncell_1st, &
54 call woverset(level, start,
end, commPressure, commVarGamma, &
55 commlamvis, commeddyvis, commpatternoverset, internaloverset)
59 do ll = 1, ntimeintervalsspectral
61 call setpointers(nn, level, ll)
63 commlamvis, commeddyvis)
72 bothpande:
if (commpressure .and. start <= irhoe .and. &
77 correctfork = getcorrectfork()
84 domains:
do nn = 1, ndom
85 do mm = 1, flowdoms(nn, level, 1)%nBocos
86 if (flowdoms(nn, level, 1)%BCType(mm) == slidinginterface)
then
90 do ll = 1, ntimeintervalsspectral
95 call setpointers(nn, level, ll)
96 call computeetotblock(icbeg(mm), icend(mm), &
97 jcbeg(mm), jcend(mm), &
98 kcbeg(mm), kcend(mm), correctfork)
109 subroutine whalo2(level, start, end, commPressure, commGamma, &
128 integer(kind=intType),
intent(in) :: level, start, end
129 logical,
intent(in) :: commPressure, commGamma, commViscous
133 integer(kind=intType) :: nn, ll
134 integer(kind=intType) :: iBeg, iEnd, jBeg, jEnd, kBeg, kEnd
136 logical :: correctForK, commLamVis, commEddyVis, commVarGamma
141 if (
viscous .and. commviscous) commlamvis = .true.
143 commeddyvis = .false.
144 if (
eddymodel .and. commviscous) commeddyvis = .true.
148 commvargamma = .false.
150 commvargamma = .true.
154 call whalo1to1(level, start,
end, commPressure, commVarGamma, &
155 commlamvis, commeddyvis, commpatterncell_2nd, &
159 call woverset(level, start,
end, commPressure, commVarGamma, &
160 commlamvis, commeddyvis, commpatternoverset, internaloverset)
164 do ll = 1, ntimeintervalsspectral
166 call setpointers(nn, level, ll)
168 commlamvis, commeddyvis)
177 bothpande:
if (commpressure .and. start <= irhoe .and. &
183 correctfork = getcorrectfork()
185 domains:
do nn = 1, ndom
190 do ll = 1, ntimeintervalsspectral
191 call setpointers(nn, level, ll)
192 call computeetotblock(2, il, 2, jl, 2, kl, correctfork)
202 calcLamVis, calcEddyVis)
217 integer(kind=intType),
intent(in) :: wstart, wend
219 logical,
intent(in) :: calcPressure, calcGamma, calcLamVis
220 logical,
intent(in) :: calcEddyVis
224 integer(kind=intType) :: oi, oj, ok, ni, nj, nk, i, l, m, n, nAvg
226 integer(kind=intType),
dimension(3) :: del
228 real(kind=realtype) :: navgreal
251 w(oi, oj, ok, l) =
zero
254 if (calcpressure)
p(oi, oj, ok) =
zero
256 if (calclamvis)
rlv(oi, oj, ok) =
zero
257 if (calceddyvis)
rev(oi, oj, ok) =
zero
263 directionloop:
do m = 1, 3
264 plusminusloop:
do i = -1, 1, 2
276 if (ni < 0 .or. ni >
ib .or. &
277 nj < 0 .or. nj >
jb .or. &
278 nk < 0 .or. nk >
kb) cycle
284 if (
iblank(ni, nj, nk) == 1)
then
292 w(oi, oj, ok, l) =
w(oi, oj, ok, l) +
w(ni, nj, nk, l)
299 p(oi, oj, ok) =
p(oi, oj, ok) +
p(ni, nj, nk)
303 rlv(oi, oj, ok) =
rlv(oi, oj, ok) +
rlv(ni, nj, nk)
305 rev(oi, oj, ok) =
rev(oi, oj, ok) +
rev(ni, nj, nk)
315 checknoneighbors:
if (navg > 0)
then
320 navgreal = real(navg, realtype)
325 w(oi, oj, ok, l) =
w(oi, oj, ok, l) / navgreal
331 if (calcpressure)
p(oi, oj, ok) =
p(oi, oj, ok) / navgreal
332 if (calcgamma)
gamma(oi, oj, ok) =
gamma(oi, oj, ok) / navgreal
333 if (calclamvis)
rlv(oi, oj, ok) =
rlv(oi, oj, ok) / navgreal
334 if (calceddyvis)
rev(oi, oj, ok) =
rev(oi, oj, ok) / navgreal
336 else checknoneighbors
342 w(oi, oj, ok, l) =
winf(l)
345 if (calcpressure)
p(oi, oj, ok) =
pinfcorr
347 if (calclamvis)
rlv(oi, oj, ok) =
muinf
350 end if checknoneighbors
357 commEddyVis, level, sps, derivPointers, nVar, varOffset)
367 integer(kind=intType),
intent(in) :: start, end, level, sps
368 logical,
intent(in) :: commPressure, commVarGamma, commLamVis, commEddyVis
369 logical,
intent(in) :: derivPointers
370 integer(kind=intType),
intent(in) :: varOffset
373 integer(kind=intType),
intent(out) :: nVar
376 integer(kind=intType) :: nn, k
377 type(
blocktype),
pointer :: blk, blk1, blkLevel
380 domainloop:
do nn = 1,
ndom
384 if (derivpointers)
then
388 blklevel =>
flowdoms(nn, level, sps)
394 blk%realCommVars(nvar)%var => blklevel%w(:, :, :, k)
397 if (commpressure)
then
399 blk%realCommVars(nvar)%var => blklevel%P(:, :, :)
402 if (commvargamma)
then
404 blk%realCommVars(nvar)%var => blk1%gamma(:, :, :)
409 blk%realCommvars(nvar)%var => blk1%rlv(:, :, :)
412 if (commeddyvis)
then
414 blk%realCommVars(nvar)%var => blklevel%rev(:, :, :)
418 nvar = nvar - varoffset
422 commVarGamma, commLamVis, commEddyVis, &
423 commPattern, internal)
441 integer(kind=intType),
intent(in) :: level, start, end
442 logical,
intent(in) :: commPressure, commVarGamma
443 logical,
intent(in) :: commLamVis, commEddyVis
445 type(
commtype),
dimension(*),
intent(in) :: commPattern
448 integer(kind=intType) :: nVar, nn, k, sps
450 logical :: correctPeriodic
456 correctperiodic = .false.
457 if (start <=
ivx .and.
end >= ivz) correctPeriodic = .true.
459 spectralmodes:
do sps = 1, ntimeintervalsspectral
462 commlamvis, commeddyvis, level, sps, .false., nvar, 0)
471 if (correctperiodic)
then
472 if (internal(level)%nPeriodic > 0)
then
474 internal(level)%nPeriodic, internal(level)%periodicData)
477 if (commpattern(level)%nPeriodic > 0)
then
479 commpattern(level)%nPeriodic, commpattern(level)%periodicData)
501 integer(kind=intType),
intent(in) :: level, sp, nPeriodic
506 integer(kind=intType) :: nn, mm, ii, i, j, k
507 real(kind=realtype) :: vx, vy, vz
509 real(kind=realtype),
dimension(3, 3) :: rotmatrix
517 rotmatrix = periodicdata(nn)%rotMatrix
521 do ii = 1, periodicdata(nn)%nhalos
525 mm = periodicdata(nn)%block(ii)
526 i = periodicdata(nn)%indices(ii, 1)
527 j = periodicdata(nn)%indices(ii, 2)
528 k = periodicdata(nn)%indices(ii, 3)
538 flowdoms(mm, level, sp)%w(i, j, k,
ivx) = rotmatrix(1, 1) * vx &
539 + rotmatrix(1, 2) * vy &
540 + rotmatrix(1, 3) * vz
541 flowdoms(mm, level, sp)%w(i, j, k,
ivy) = rotmatrix(2, 1) * vx &
542 + rotmatrix(2, 2) * vy &
543 + rotmatrix(2, 3) * vz
544 flowdoms(mm, level, sp)%w(i, j, k,
ivz) = rotmatrix(3, 1) * vx &
545 + rotmatrix(3, 2) * vy &
546 + rotmatrix(3, 3) * vz
568 integer(kind=intType),
intent(in) :: level, sps
570 type(
commtype),
dimension(*),
intent(in) :: commPattern
576 integer :: size, procID, ierr, index
577 integer,
dimension(mpi_status_size) :: mpiStatus
579 integer(kind=intType) :: nVar, mm
580 integer(kind=intType) :: i, j, k, ii, jj
581 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
587 sends:
do i = 1, commpattern(level)%nProcSend
592 procid = commpattern(level)%sendProc(i)
593 size = nvar * commpattern(level)%nsend(i)
599 do j = 1, commpattern(level)%nsend(i)
604 d1 = commpattern(level)%sendList(i)%block(j)
605 i1 = commpattern(level)%sendList(i)%indices(j, 1) + 1
606 j1 = commpattern(level)%sendList(i)%indices(j, 2) + 1
607 k1 = commpattern(level)%sendList(i)%indices(j, 3) + 1
620 call mpi_isend(
sendbuffer(ii),
size, adflow_real, procid, &
623 call echk(ierr, __file__, __line__)
634 receives:
do i = 1, commpattern(level)%nProcRecv
639 procid = commpattern(level)%recvProc(i)
640 size = nvar * commpattern(level)%nrecv(i)
644 call mpi_irecv(
recvbuffer(ii),
size, adflow_real, procid, &
646 call echk(ierr, __file__, __line__)
657 localcopy:
do i = 1, internal(level)%ncopy
661 d1 = internal(level)%donorBlock(i)
662 i1 = internal(level)%donorIndices(i, 1) + 1
663 j1 = internal(level)%donorIndices(i, 2) + 1
664 k1 = internal(level)%donorIndices(i, 3) + 1
668 d2 = internal(level)%haloBlock(i)
669 i2 = internal(level)%haloIndices(i, 1) + 1
670 j2 = internal(level)%haloIndices(i, 2) + 1
671 k2 = internal(level)%haloIndices(i, 3) + 1
674 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) = &
675 flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1)
683 size = commpattern(level)%nProcRecv
684 completerecvs:
do i = 1, commpattern(level)%nProcRecv
688 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
689 call echk(ierr, __file__, __line__)
694 jj = nvar * commpattern(level)%nrecvCum(ii - 1)
696 do j = 1, commpattern(level)%nrecv(ii)
700 d2 = commpattern(level)%recvList(ii)%block(j)
701 i2 = commpattern(level)%recvList(ii)%indices(j, 1) + 1
702 j2 = commpattern(level)%recvList(ii)%indices(j, 2) + 1
703 k2 = commpattern(level)%recvList(ii)%indices(j, 3) + 1
714 size = commpattern(level)%nProcSend
715 do i = 1, commpattern(level)%nProcSend
716 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
735 integer(kind=intType),
intent(in) :: level, sps
737 type(
commtype),
dimension(*),
intent(in) :: commPattern
743 integer :: size, procID, ierr, index
744 integer,
dimension(mpi_status_size) :: mpiStatus
746 integer(kind=intType) :: nVar, mm
747 integer(kind=intType) :: i, j, k, ii, jj
748 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
756 recvs:
do i = 1, commpattern(level)%nProcRecv
761 procid = commpattern(level)%recvProc(i)
762 size = nvar * commpattern(level)%nrecv(i)
766 do j = 1, commpattern(level)%nrecv(i)
770 d2 = commpattern(level)%recvList(i)%block(j)
771 i2 = commpattern(level)%recvList(i)%indices(j, 1) + 1
772 j2 = commpattern(level)%recvList(i)%indices(j, 2) + 1
773 k2 = commpattern(level)%recvList(i)%indices(j, 3) + 1
777 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) =
zero
783 call mpi_isend(
recvbuffer(ii),
size, adflow_real, procid, &
786 call echk(ierr, __file__, __line__)
797 sends:
do i = 1, commpattern(level)%nProcSend
802 procid = commpattern(level)%sendProc(i)
803 size = nvar * commpattern(level)%nsend(i)
807 call mpi_irecv(
sendbuffer(ii),
size, adflow_real, procid, &
809 call echk(ierr, __file__, __line__)
819 localcopy:
do i = 1, internal(level)%ncopy
823 d1 = internal(level)%donorBlock(i)
824 i1 = internal(level)%donorIndices(i, 1) + 1
825 j1 = internal(level)%donorIndices(i, 2) + 1
826 k1 = internal(level)%donorIndices(i, 3) + 1
830 d2 = internal(level)%haloBlock(i)
831 i2 = internal(level)%haloIndices(i, 1) + 1
832 j2 = internal(level)%haloIndices(i, 2) + 1
833 k2 = internal(level)%haloIndices(i, 3) + 1
837 flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) = &
838 flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) + &
839 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2)
840 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) =
zero
848 size = commpattern(level)%nProcSend
849 completesends:
do i = 1, commpattern(level)%nProcSend
853 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
854 call echk(ierr, __file__, __line__)
859 jj = nvar * commpattern(level)%nsendCum(ii - 1)
861 do j = 1, commpattern(level)%nsend(ii)
865 d2 = commpattern(level)%sendList(ii)%block(j)
866 i2 = commpattern(level)%sendList(ii)%indices(j, 1) + 1
867 j2 = commpattern(level)%sendList(ii)%indices(j, 2) + 1
868 k2 = commpattern(level)%sendList(ii)%indices(j, 3) + 1
874 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) = &
883 size = commpattern(level)%nProcRecv
884 do i = 1, commpattern(level)%nProcRecv
885 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
886 call echk(ierr, __file__, __line__)
906 integer(kind=intType),
intent(in) :: level, sps
908 type(
commtype),
dimension(*),
intent(in) :: commPattern
914 integer :: size, procID, ierr, index
915 integer,
dimension(mpi_status_size) :: mpiStatus
917 integer(kind=intType) :: nVar, mm
918 integer(kind=intType) :: i, j, k, ii, jj
919 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
921 integer(kind=intType),
dimension(:),
allocatable :: sendBufInt
922 integer(kind=intType),
dimension(:),
allocatable :: recvBufInt
924 ii = commpattern(level)%nProcSend
925 ii = commpattern(level)%nsendCum(ii)
926 jj = commpattern(level)%nProcRecv
927 jj = commpattern(level)%nrecvCum(jj)
929 allocate (sendbufint(ii * nvar), recvbufint(jj * nvar), stat=ierr)
935 sends:
do i = 1, commpattern(level)%nProcSend
940 procid = commpattern(level)%sendProc(i)
941 size = nvar * commpattern(level)%nsend(i)
947 do j = 1, commpattern(level)%nsend(i)
952 d1 = commpattern(level)%sendList(i)%block(j)
953 i1 = commpattern(level)%sendList(i)%indices(j, 1) + 1
954 j1 = commpattern(level)%sendList(i)%indices(j, 2) + 1
955 k1 = commpattern(level)%sendList(i)%indices(j, 3) + 1
961 sendbufint(jj) =
flowdoms(d1, level, sps)%intCommVars(k)%var(i1, j1, k1)
968 call mpi_isend(sendbufint(ii),
size, adflow_integer, procid, &
971 call echk(ierr, __file__, __line__)
982 receives:
do i = 1, commpattern(level)%nProcRecv
987 procid = commpattern(level)%recvProc(i)
988 size = nvar * commpattern(level)%nrecv(i)
992 call mpi_irecv(recvbufint(ii),
size, adflow_integer, procid, &
994 call echk(ierr, __file__, __line__)
1004 localcopy:
do i = 1, internal(level)%ncopy
1008 d1 = internal(level)%donorBlock(i)
1009 i1 = internal(level)%donorIndices(i, 1) + 1
1010 j1 = internal(level)%donorIndices(i, 2) + 1
1011 k1 = internal(level)%donorIndices(i, 3) + 1
1015 d2 = internal(level)%haloBlock(i)
1016 i2 = internal(level)%haloIndices(i, 1) + 1
1017 j2 = internal(level)%haloIndices(i, 2) + 1
1018 k2 = internal(level)%haloIndices(i, 3) + 1
1021 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) = &
1022 flowdoms(d1, level, sps)%intCommVars(k)%var(i1, j1, k1)
1030 size = commpattern(level)%nProcRecv
1031 completerecvs:
do i = 1, commpattern(level)%nProcRecv
1035 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
1036 call echk(ierr, __file__, __line__)
1041 jj = nvar * commpattern(level)%nrecvCum(ii - 1)
1043 do j = 1, commpattern(level)%nrecv(ii)
1047 d2 = commpattern(level)%recvList(ii)%block(j)
1048 i2 = commpattern(level)%recvList(ii)%indices(j, 1) + 1
1049 j2 = commpattern(level)%recvList(ii)%indices(j, 2) + 1
1050 k2 = commpattern(level)%recvList(ii)%indices(j, 3) + 1
1054 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) = recvbufint(jj)
1057 end do completerecvs
1061 size = commpattern(level)%nProcSend
1062 do i = 1, commpattern(level)%nProcSend
1063 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
1064 call echk(ierr, __file__, __line__)
1067 deallocate (recvbufint, sendbufint)
1086 integer(kind=intType),
intent(in) :: level, sps
1088 type(
commtype),
dimension(*),
intent(in) :: commPattern
1094 integer :: size, procID, ierr, index
1095 integer,
dimension(mpi_status_size) :: mpiStatus
1097 integer(kind=intType) :: nVar, mm
1098 integer(kind=intType) :: i, j, k, ii, jj
1099 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
1100 integer(kind=intType),
dimension(:),
allocatable :: sendBufInt
1101 integer(kind=intType),
dimension(:),
allocatable :: recvBufInt
1103 ii = commpattern(level)%nProcSend
1104 ii = commpattern(level)%nsendCum(ii)
1105 jj = commpattern(level)%nProcRecv
1106 jj = commpattern(level)%nrecvCum(jj)
1108 allocate (sendbufint(ii * nvar), recvbufint(jj * nvar), stat=ierr)
1116 recvs:
do i = 1, commpattern(level)%nProcRecv
1121 procid = commpattern(level)%recvProc(i)
1122 size = nvar * commpattern(level)%nrecv(i)
1126 do j = 1, commpattern(level)%nrecv(i)
1130 d2 = commpattern(level)%recvList(i)%block(j)
1131 i2 = commpattern(level)%recvList(i)%indices(j, 1) + 1
1132 j2 = commpattern(level)%recvList(i)%indices(j, 2) + 1
1133 k2 = commpattern(level)%recvList(i)%indices(j, 3) + 1
1136 recvbufint(jj) =
flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2)
1137 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) = 0
1143 call mpi_isend(recvbufint(ii),
size, adflow_integer, procid, &
1146 call echk(ierr, __file__, __line__)
1157 sends:
do i = 1, commpattern(level)%nProcSend
1162 procid = commpattern(level)%sendProc(i)
1163 size = nvar * commpattern(level)%nsend(i)
1167 call mpi_irecv(sendbufint(ii),
size, adflow_integer, procid, &
1169 call echk(ierr, __file__, __line__)
1179 localcopy:
do i = 1, internal(level)%ncopy
1183 d1 = internal(level)%donorBlock(i)
1184 i1 = internal(level)%donorIndices(i, 1) + 1
1185 j1 = internal(level)%donorIndices(i, 2) + 1
1186 k1 = internal(level)%donorIndices(i, 3) + 1
1190 d2 = internal(level)%haloBlock(i)
1191 i2 = internal(level)%haloIndices(i, 1) + 1
1192 j2 = internal(level)%haloIndices(i, 2) + 1
1193 k2 = internal(level)%haloIndices(i, 3) + 1
1198 flowdoms(d1, level, sps)%intCommVars(k)%var(i1, j1, k1) = &
1199 flowdoms(d1, level, sps)%intCommVars(k)%var(i1, j1, k1) + &
1200 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2)
1201 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) = 0
1209 size = commpattern(level)%nProcSend
1210 completesends:
do i = 1, commpattern(level)%nProcSend
1214 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
1215 call echk(ierr, __file__, __line__)
1220 jj = nvar * commpattern(level)%nsendCum(ii - 1)
1222 do j = 1, commpattern(level)%nsend(ii)
1226 d2 = commpattern(level)%sendList(ii)%block(j)
1227 i2 = commpattern(level)%sendList(ii)%indices(j, 1) + 1
1228 j2 = commpattern(level)%sendList(ii)%indices(j, 2) + 1
1229 k2 = commpattern(level)%sendList(ii)%indices(j, 3) + 1
1235 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) = &
1236 flowdoms(d2, level, sps)%intCommVars(k)%var(i2, j2, k2) + sendbufint(jj)
1240 end do completesends
1244 size = commpattern(level)%nProcRecv
1245 do i = 1, commpattern(level)%nProcRecv
1246 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
1247 call echk(ierr, __file__, __line__)
1250 deallocate (recvbufint, sendbufint)
1255 commVarGamma, commLamVis, commEddyVis, &
1256 commPattern, internal)
1267 integer(kind=intType),
intent(in) :: level, start, end
1268 logical,
intent(in) :: commPressure, commVarGamma
1269 logical,
intent(in) :: commLamVis, commEddyVis
1271 type(
commtype),
dimension(*),
intent(in) :: commPattern
1274 integer(kind=intType) :: nVar, nn, k, sps
1279 commlamvis, commeddyvis, level, sps, .true., nvar, 0)
1288 end do spectralmodes
1293 commVarGamma, commLamVis, commEddyVis, &
1294 commPattern, internal)
1305 integer(kind=intType),
intent(in) :: level, start, end
1306 logical,
intent(in) :: commPressure, commVarGamma
1307 logical,
intent(in) :: commLamVis, commEddyVis
1309 type(
commtype),
dimension(*),
intent(in) :: commPattern
1312 integer(kind=intType) :: nVar, nn, k, sps
1317 commlamvis, commeddyvis, level, sps, .true., nvar, 0)
1326 end do spectralmodes
1331 commVarGamma, commLamVis, commEddyVis, &
1332 commPattern, internal)
1351 integer(kind=intType),
intent(in) :: level, start, end
1352 logical,
intent(in) :: commPressure, commVarGamma
1353 logical,
intent(in) :: commLamVis, commEddyVis
1355 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1359 integer(kind=intType) :: nVar, sps
1364 commlamvis, commeddyvis, level, sps, .false., nvar, 0)
1372 end do spectralmodes
1377 commVarGamma, commLamVis, commEddyVis, &
1378 commPattern, internal)
1397 integer(kind=intType),
intent(in) :: level, start, end
1398 logical,
intent(in) :: commPressure, commVarGamma
1399 logical,
intent(in) :: commLamVis, commEddyVis
1401 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1403 integer(kind=intType) :: nVar, sps, offset
1410 commlamvis, commeddyvis, level, sps, .true., nvar, 0)
1415 commlamvis, commeddyvis, level, sps, .false., nvar, offset)
1423 end do spectralmodes
1427 commVarGamma, commLamVis, commEddyVis, &
1428 commPattern, internal)
1441 integer(kind=intType),
intent(in) :: level, start, end
1442 logical,
intent(in) :: commPressure, commVarGamma
1443 logical,
intent(in) :: commLamVis, commEddyVis
1445 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1447 integer(kind=intType) :: nVar, sps, offset
1454 commlamvis, commeddyvis, level, sps, .true., nvar, 0)
1459 commlamvis, commeddyvis, level, sps, .false., nvar, offset)
1467 end do spectralmodes
1484 integer(kind=intType),
intent(in) :: level, sps
1486 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1491 integer :: size, procId, ierr, index
1492 integer,
dimension(mpi_status_size) :: mpiStatus
1494 integer(kind=intType) :: nVar
1495 integer(kind=intType) :: i, j, k, ii, jj
1496 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
1497 real(kind=realtype),
dimension(:),
pointer :: weight
1503 sends:
do i = 1, commpattern(level, sps)%nProcSend
1508 procid = commpattern(level, sps)%sendProc(i)
1509 size = nvar * commpattern(level, sps)%nsend(i)
1515 do j = 1, commpattern(level, sps)%nsend(i)
1520 d1 = commpattern(level, sps)%sendList(i)%block(j)
1521 i1 = commpattern(level, sps)%sendList(i)%indices(j, 1) + 1
1522 j1 = commpattern(level, sps)%sendList(i)%indices(j, 2) + 1
1523 k1 = commpattern(level, sps)%sendList(i)%indices(j, 3) + 1
1524 weight => commpattern(level, sps)%sendList(i)%interp(j, :)
1531 weight(1) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) + &
1532 weight(2) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1) + &
1533 weight(3) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1) + &
1534 weight(4) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1) + &
1535 weight(5) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1 + 1) + &
1536 weight(6) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1 + 1) + &
1537 weight(7) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1 + 1) + &
1538 weight(8) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1 + 1)
1545 call mpi_isend(
sendbuffer(ii),
size, adflow_real, procid, &
1548 call echk(ierr, __file__, __line__)
1559 receives:
do i = 1, commpattern(level, sps)%nProcRecv
1564 procid = commpattern(level, sps)%recvProc(i)
1565 size = nvar * commpattern(level, sps)%nrecv(i)
1569 call mpi_irecv(
recvbuffer(ii),
size, adflow_real, procid, &
1571 call echk(ierr, __file__, __line__)
1581 localinterp:
do i = 1, internal(level, sps)%ncopy
1585 d1 = internal(level, sps)%donorBlock(i)
1586 i1 = internal(level, sps)%donorIndices(i, 1) + 1
1587 j1 = internal(level, sps)%donorIndices(i, 2) + 1
1588 k1 = internal(level, sps)%donorIndices(i, 3) + 1
1590 weight => internal(level, sps)%donorInterp(i, :)
1594 d2 = internal(level, sps)%haloBlock(i)
1595 i2 = internal(level, sps)%haloIndices(i, 1) + 1
1596 j2 = internal(level, sps)%haloIndices(i, 2) + 1
1597 k2 = internal(level, sps)%haloIndices(i, 3) + 1
1602 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) = &
1603 weight(1) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) + &
1604 weight(2) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1) + &
1605 weight(3) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1) + &
1606 weight(4) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1) + &
1607 weight(5) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1 + 1) + &
1608 weight(6) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1 + 1) + &
1609 weight(7) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1 + 1) + &
1610 weight(8) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1 + 1)
1617 size = commpattern(level, sps)%nProcRecv
1618 completerecvs:
do i = 1, commpattern(level, sps)%nProcRecv
1622 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
1623 call echk(ierr, __file__, __line__)
1628 jj = nvar * commpattern(level, sps)%nrecvCum(ii - 1)
1630 do j = 1, commpattern(level, sps)%nrecv(ii)
1634 d2 = commpattern(level, sps)%recvList(ii)%block(j)
1635 i2 = commpattern(level, sps)%recvList(ii)%indices(j, 1) + 1
1636 j2 = commpattern(level, sps)%recvList(ii)%indices(j, 2) + 1
1637 k2 = commpattern(level, sps)%recvList(ii)%indices(j, 3) + 1
1644 end do completerecvs
1648 size = commpattern(level, sps)%nProcSend
1649 do i = 1, commpattern(level, sps)%nProcSend
1650 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
1651 call echk(ierr, __file__, __line__)
1669 integer(kind=intType),
intent(in) :: level, sps
1671 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1676 integer :: size, procId, ierr, index
1677 integer,
dimension(mpi_status_size) :: mpiStatus
1679 integer(kind=intType) :: nVar
1680 integer(kind=intType) :: i, j, k, ii, jj
1681 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
1682 real(kind=realtype),
dimension(:),
pointer :: weight, weightd
1688 sends:
do i = 1, commpattern(level, sps)%nProcSend
1693 procid = commpattern(level, sps)%sendProc(i)
1694 size = nvar * commpattern(level, sps)%nsend(i)
1700 do j = 1, commpattern(level, sps)%nsend(i)
1705 d1 = commpattern(level, sps)%sendList(i)%block(j)
1706 i1 = commpattern(level, sps)%sendList(i)%indices(j, 1) + 1
1707 j1 = commpattern(level, sps)%sendList(i)%indices(j, 2) + 1
1708 k1 = commpattern(level, sps)%sendList(i)%indices(j, 3) + 1
1709 weight => commpattern(level, sps)%sendList(i)%interp(j, :)
1710 weightd => commpattern(level, sps)%sendList(i)%interpd(j, :)
1717 weight(1) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) + &
1718 weight(2) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1) + &
1719 weight(3) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1) + &
1720 weight(4) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1) + &
1721 weight(5) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1 + 1) + &
1722 weight(6) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1 + 1) + &
1723 weight(7) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1 + 1) + &
1724 weight(8) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1 + 1) + &
1725 weightd(1) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1, k1) + &
1726 weightd(2) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1, k1) + &
1727 weightd(3) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1 + 1, k1) + &
1728 weightd(4) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1 + 1, k1) + &
1729 weightd(5) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1, k1 + 1) + &
1730 weightd(6) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1, k1 + 1) + &
1731 weightd(7) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1 + 1, k1 + 1) + &
1732 weightd(8) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1 + 1, k1 + 1)
1740 call mpi_isend(
sendbuffer(ii),
size, adflow_real, procid, &
1743 call echk(ierr, __file__, __line__)
1754 receives:
do i = 1, commpattern(level, sps)%nProcRecv
1759 procid = commpattern(level, sps)%recvProc(i)
1760 size = nvar * commpattern(level, sps)%nrecv(i)
1764 call mpi_irecv(
recvbuffer(ii),
size, adflow_real, procid, &
1766 call echk(ierr, __file__, __line__)
1776 localinterp:
do i = 1, internal(level, sps)%ncopy
1780 d1 = internal(level, sps)%donorBlock(i)
1781 i1 = internal(level, sps)%donorIndices(i, 1) + 1
1782 j1 = internal(level, sps)%donorIndices(i, 2) + 1
1783 k1 = internal(level, sps)%donorIndices(i, 3) + 1
1785 weight => internal(level, sps)%donorInterp(i, :)
1786 weightd => internal(level, sps)%donorInterpd(i, :)
1790 d2 = internal(level, sps)%haloBlock(i)
1791 i2 = internal(level, sps)%haloIndices(i, 1) + 1
1792 j2 = internal(level, sps)%haloIndices(i, 2) + 1
1793 k2 = internal(level, sps)%haloIndices(i, 3) + 1
1798 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) = &
1799 weight(1) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1) + &
1800 weight(2) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1) + &
1801 weight(3) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1) + &
1802 weight(4) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1) + &
1803 weight(5) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1, k1 + 1) + &
1804 weight(6) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1, k1 + 1) + &
1805 weight(7) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1, j1 + 1, k1 + 1) + &
1806 weight(8) *
flowdoms(d1, level, sps)%realCommVars(k)%var(i1 + 1, j1 + 1, k1 + 1) + &
1807 weightd(1) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1, k1) + &
1808 weightd(2) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1, k1) + &
1809 weightd(3) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1 + 1, k1) + &
1810 weightd(4) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1 + 1, k1) + &
1811 weightd(5) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1, k1 + 1) + &
1812 weightd(6) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1, k1 + 1) + &
1813 weightd(7) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1, j1 + 1, k1 + 1) + &
1814 weightd(8) *
flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(i1 + 1, j1 + 1, k1 + 1)
1822 size = commpattern(level, sps)%nProcRecv
1823 completerecvs:
do i = 1, commpattern(level, sps)%nProcRecv
1827 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
1828 call echk(ierr, __file__, __line__)
1833 jj = nvar * commpattern(level, sps)%nrecvCum(ii - 1)
1835 do j = 1, commpattern(level, sps)%nrecv(ii)
1839 d2 = commpattern(level, sps)%recvList(ii)%block(j)
1840 i2 = commpattern(level, sps)%recvList(ii)%indices(j, 1) + 1
1841 j2 = commpattern(level, sps)%recvList(ii)%indices(j, 2) + 1
1842 k2 = commpattern(level, sps)%recvList(ii)%indices(j, 3) + 1
1849 end do completerecvs
1853 size = commpattern(level, sps)%nProcSend
1854 do i = 1, commpattern(level, sps)%nProcSend
1855 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
1856 call echk(ierr, __file__, __line__)
1874 integer(kind=intType),
intent(in) :: level, sps
1876 type(
commtype),
dimension(:, :),
intent(in) :: commPattern
1881 integer :: size, procId, ierr, index
1882 integer,
dimension(mpi_status_size) :: mpiStatus
1884 integer(kind=intType) :: nVar
1885 integer(kind=intType) :: i, j, k, ii, jj, kk
1886 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2, iii, jjj, kkk
1887 real(kind=realtype),
dimension(:),
pointer :: weight, weightd
1888 real(kind=realtype) :: vard
1894 recvs:
do i = 1, commpattern(level, sps)%nProcRecv
1899 procid = commpattern(level, sps)%recvProc(i)
1900 size = nvar * commpattern(level, sps)%nrecv(i)
1904 do j = 1, commpattern(level, sps)%nrecv(i)
1908 d2 = commpattern(level, sps)%recvList(i)%block(j)
1909 i2 = commpattern(level, sps)%recvList(i)%indices(j, 1) + 1
1910 j2 = commpattern(level, sps)%recvList(i)%indices(j, 2) + 1
1911 k2 = commpattern(level, sps)%recvList(i)%indices(j, 3) + 1
1916 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) =
zero
1921 call mpi_isend(
recvbuffer(ii),
size, adflow_real, procid, &
1924 call echk(ierr, __file__, __line__)
1935 sends:
do i = 1, commpattern(level, sps)%nProcSend
1940 procid = commpattern(level, sps)%sendProc(i)
1941 size = nvar * commpattern(level, sps)%nsend(i)
1945 call mpi_irecv(
sendbuffer(ii),
size, adflow_real, procid, &
1947 call echk(ierr, __file__, __line__)
1957 localinterp:
do i = 1, internal(level, sps)%ncopy
1961 d1 = internal(level, sps)%donorBlock(i)
1962 i1 = internal(level, sps)%donorIndices(i, 1) + 1
1963 j1 = internal(level, sps)%donorIndices(i, 2) + 1
1964 k1 = internal(level, sps)%donorIndices(i, 3) + 1
1966 weight => internal(level, sps)%donorInterp(i, :)
1967 weightd => internal(level, sps)%donorInterpd(i, :)
1971 d2 = internal(level, sps)%haloBlock(i)
1972 i2 = internal(level, sps)%haloIndices(i, 1) + 1
1973 j2 = internal(level, sps)%haloIndices(i, 2) + 1
1974 k2 = internal(level, sps)%haloIndices(i, 3) + 1
1979 vard =
flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2)
1985 flowdoms(d1, level, sps)%realCommVars(k)%var(iii, jjj, kkk) = &
1986 flowdoms(d1, level, sps)%realCommVars(k)%var(iii, jjj, kkk) + &
1989 weightd(kk) = weightd(kk) + &
1990 flowdoms(d1, level, sps)%realCommVars(k + nvar)%var(iii, jjj, kkk) * vard
1995 flowdoms(d2, level, sps)%realCommVars(k)%var(i2, j2, k2) =
zero
2002 size = commpattern(level, sps)%nProcSend
2003 completesends:
do i = 1, commpattern(level, sps)%nProcSend
2007 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
2008 call echk(ierr, __file__, __line__)
2014 jj = nvar * commpattern(level, sps)%nsendCum(ii - 1)
2016 do j = 1, commpattern(level, sps)%nsend(ii)
2020 d2 = commpattern(level, sps)%sendList(ii)%block(j)
2021 i2 = commpattern(level, sps)%sendList(ii)%indices(j, 1) + 1
2022 j2 = commpattern(level, sps)%sendList(ii)%indices(j, 2) + 1
2023 k2 = commpattern(level, sps)%sendList(ii)%indices(j, 3) + 1
2025 weight => commpattern(level, sps)%sendList(ii)%interp(j, :)
2026 weightd => commpattern(level, sps)%sendList(ii)%interpd(j, :)
2037 flowdoms(d2, level, sps)%realCommVars(k)%var(iii, jjj, kkk) = &
2038 flowdoms(d2, level, sps)%realCommVars(k)%var(iii, jjj, kkk) + &
2041 weightd(kk) = weightd(kk) + &
2042 flowdoms(d2, level, sps)%realCommVars(k + nvar)%var(iii, jjj, kkk) * vard
2050 end do completesends
2054 size = commpattern(level, sps)%nProcRecv
2055 do i = 1, commpattern(level, sps)%nProcRecv
2056 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
2057 call echk(ierr, __file__, __line__)
2063 subroutine whalo2_b(level, start, end, commPressure, commGamma, &
2082 integer(kind=intType),
intent(in) :: level, start, end
2083 logical,
intent(in) :: commPressure, commGamma, commViscous
2087 integer(kind=intType) :: nn, mm, ll
2088 integer(kind=intType) :: iBeg, iEnd, jBeg, jEnd, kBeg, kEnd
2090 logical :: correctForK, commLamVis, commEddyVis, commVarGamma
2094 commlamvis = .false.
2095 if (
viscous .and. commviscous) commlamvis = .true.
2097 commeddyvis = .false.
2098 if (
eddymodel .and. commviscous) commeddyvis = .true.
2102 commvargamma = .false.
2104 commvargamma = .true.
2106 bothpande:
if (commpressure .and. start <=
irhoe .and. &
2112 correctfork = getcorrectfork()
2114 domains:
do nn = 1, ndom
2119 do ll = 1, ntimeintervalsspectral
2120 call setpointers_b(nn, level, ll)
2121 call computeetotblock_b(2, il, 2, jl, 2, kl, correctfork)
2126 call woverset_b(level, start,
end, commPressure, commVarGamma, &
2127 commlamvis, commeddyvis, commpatternoverset, internaloverset)
2130 call whalo1to1_b(level, start,
end, commPressure, commVarGamma, &
2131 commlamvis, commeddyvis, commpatterncell_2nd, &
2140 subroutine whalo2_d(level, start, end, commPressure, commGamma, &
2159 integer(kind=intType),
intent(in) :: level, start, end
2160 logical,
intent(in) :: commPressure, commGamma, commViscous
2164 integer(kind=intType) :: nn, mm, ll
2165 integer(kind=intType) :: iBeg, iEnd, jBeg, jEnd, kBeg, kEnd
2167 logical :: correctForK, commLamVis, commEddyVis, commVarGamma
2171 commlamvis = .false.
2172 if (
viscous .and. commviscous) commlamvis = .true.
2174 commeddyvis = .false.
2175 if (
eddymodel .and. commviscous) commeddyvis = .true.
2179 commvargamma = .false.
2181 commvargamma = .true.
2185 call whalo1to1_d(level, start,
end, commPressure, commVarGamma, &
2186 commlamvis, commeddyvis, commpatterncell_2nd, &
2190 call woverset_d(level, start,
end, commPressure, commVarGamma, &
2191 commlamvis, commeddyvis, commpatternoverset, internaloverset)
2202 bothpande:
if (commpressure .and. start <= irhoe .and. &
2208 correctfork = getcorrectfork()
2210 domains:
do nn = 1, ndom
2215 do ll = 1, ntimeintervalsspectral
2216 call setpointers_d(nn, level, ll)
2217 call computeetotblock_d(2, il, 2, jl, 2, kl, correctfork)
2241 integer(kind=intType) :: level, start, end
2245 integer :: size, procID, ierr, index
2246 integer,
dimension(mpi_status_size) :: mpiStatus
2248 integer(kind=intType) :: nVar, sps
2249 integer(kind=intType) :: ii, jj, mm, nn, i, j, k, l
2250 integer(kind=intType) :: dd1, ii1, jj1, kk1, dd2, ii2, jj2, kk2
2252 real(kind=realtype),
pointer,
dimension(:, :, :) :: ddw1, ddw2
2256 nvar = max(0_inttype, (
end - start + 1))
2257 if (nvar == 0)
return
2262 spectralloop:
do sps = 1, ntimeintervalsspectral
2263 domains:
do mm = 1, ndom
2267 call setpointers(mm, level, sps)
2274 bocos:
do nn = 1, nbocos
2279 select case (bcfaceid(nn))
2281 ddw1 => dw(1, 1:, 1:, :); ddw2 => dw(2, 1:, 1:, :)
2283 ddw1 => dw(ie, 1:, 1:, :); ddw2 => dw(il, 1:, 1:, :)
2285 ddw1 => dw(1:, 1, 1:, :); ddw2 => dw(1:, 2, 1:, :)
2287 ddw1 => dw(1:, je, 1:, :); ddw2 => dw(1:, jl, 1:, :)
2289 ddw1 => dw(1:, 1:, 1, :); ddw2 => dw(1:, 1:, 2, :)
2291 ddw1 => dw(1:, 1:, ke, :); ddw2 => dw(1:, 1:, kl, :)
2301 ddw1(i, j, l) = ddw2(i, j, l)
2313 sends:
do i = 1, commpatterncell_1st(level)%nProcSend
2318 procid = commpatterncell_1st(level)%sendProc(i)
2319 size = nvar * commpatterncell_1st(level)%nsend(i)
2325 do j = 1, commpatterncell_1st(level)%nsend(i)
2330 dd1 = commpatterncell_1st(level)%sendList(i)%block(j)
2331 ii1 = commpatterncell_1st(level)%sendList(i)%indices(j, 1)
2332 jj1 = commpatterncell_1st(level)%sendList(i)%indices(j, 2)
2333 kk1 = commpatterncell_1st(level)%sendList(i)%indices(j, 3)
2339 sendbuffer(jj) = flowdoms(dd1, level, sps)%dw(ii1, jj1, kk1, k)
2347 call mpi_isend(sendbuffer(ii),
size, adflow_real, procid, &
2348 procid, adflow_comm_world, sendrequests(i), &
2350 call echk(ierr, __file__, __line__)
2361 receives:
do i = 1, commpatterncell_1st(level)%nProcRecv
2366 procid = commpatterncell_1st(level)%recvProc(i)
2367 size = nvar * commpatterncell_1st(level)%nrecv(i)
2371 call mpi_irecv(recvbuffer(ii),
size, adflow_real, procid, &
2372 myid, adflow_comm_world, recvrequests(i), ierr)
2373 call echk(ierr, __file__, __line__)
2383 localcopy:
do i = 1, internalcell_1st(level)%ncopy
2387 dd1 = internalcell_1st(level)%donorBlock(i)
2388 ii1 = internalcell_1st(level)%donorIndices(i, 1)
2389 jj1 = internalcell_1st(level)%donorIndices(i, 2)
2390 kk1 = internalcell_1st(level)%donorIndices(i, 3)
2394 dd2 = internalcell_1st(level)%haloBlock(i)
2395 ii2 = internalcell_1st(level)%haloIndices(i, 1)
2396 jj2 = internalcell_1st(level)%haloIndices(i, 2)
2397 kk2 = internalcell_1st(level)%haloIndices(i, 3)
2402 flowdoms(dd2, level, sps)%dw(ii2, jj2, kk2, k) = &
2403 flowdoms(dd1, level, sps)%dw(ii1, jj1, kk1, k)
2411 size = commpatterncell_1st(level)%nProcRecv
2412 completerecvs:
do i = 1, commpatterncell_1st(level)%nProcRecv
2416 call mpi_waitany(
size, recvrequests, index, mpistatus, ierr)
2417 call echk(ierr, __file__, __line__)
2422 jj = nvar * commpatterncell_1st(level)%nrecvCum(ii - 1) + 1
2424 do j = 1, commpatterncell_1st(level)%nrecv(ii)
2428 dd2 = commpatterncell_1st(level)%recvList(ii)%block(j)
2429 ii2 = commpatterncell_1st(level)%recvList(ii)%indices(j, 1)
2430 jj2 = commpatterncell_1st(level)%recvList(ii)%indices(j, 2)
2431 kk2 = commpatterncell_1st(level)%recvList(ii)%indices(j, 3)
2436 flowdoms(dd2, level, sps)%dw(ii2, jj2, kk2, k) = recvbuffer(jj)
2442 end do completerecvs
2446 size = commpatterncell_1st(level)%nProcSend
2447 do i = 1, commpatterncell_1st(level)%nProcSend
2448 call mpi_waitany(
size, sendrequests, index, mpistatus, ierr)
2449 call echk(ierr, __file__, __line__)
2468 integer(kind=intType),
intent(in) :: level
2472 integer :: size, procID, ierr, index
2473 integer,
dimension(mpi_status_size) :: mpiStatus
2475 integer(kind=intType) :: i, j, ii, jj, mm
2476 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
2520 call mpi_isend(
sendbuffer(ii),
size, adflow_real, procid, &
2523 call echk(ierr, __file__, __line__)
2544 call mpi_irecv(
recvbuffer(ii),
size, adflow_real, procid, &
2546 call echk(ierr, __file__, __line__)
2571 flowdoms(d2, level, mm)%x(i2, j2, k2, 1) = &
2572 flowdoms(d1, level, mm)%x(i1, j1, k1, 1)
2573 flowdoms(d2, level, mm)%x(i2, j2, k2, 2) = &
2574 flowdoms(d1, level, mm)%x(i1, j1, k1, 2)
2575 flowdoms(d2, level, mm)%x(i2, j2, k2, 3) = &
2576 flowdoms(d1, level, mm)%x(i1, j1, k1, 3)
2595 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
2596 call echk(ierr, __file__, __line__)
2621 end do completerecvs
2634 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
2635 call echk(ierr, __file__, __line__)
2655 integer(kind=intType),
intent(in) :: level, sp, nPeriodic
2660 integer(kind=intType) :: nn, mm, ii, i, j, k
2661 real(kind=realtype) :: dx, dy, dz
2663 real(kind=realtype),
dimension(3, 3) :: rotmatrix
2664 real(kind=realtype),
dimension(3) :: rotcenter, translation
2668 do nn = 1, nperiodic
2673 rotmatrix = periodicdata(nn)%rotMatrix
2674 rotcenter = periodicdata(nn)%rotCenter
2675 translation = periodicdata(nn)%translation + rotcenter
2679 do ii = 1, periodicdata(nn)%nHalos
2683 mm = periodicdata(nn)%block(ii)
2684 i = periodicdata(nn)%indices(ii, 1)
2685 j = periodicdata(nn)%indices(ii, 2)
2686 k = periodicdata(nn)%indices(ii, 3)
2691 dx =
flowdoms(mm, level, sp)%x(i, j, k, 1) - rotcenter(1)
2692 dy =
flowdoms(mm, level, sp)%x(i, j, k, 2) - rotcenter(2)
2693 dz =
flowdoms(mm, level, sp)%x(i, j, k, 3) - rotcenter(3)
2697 flowdoms(mm, level, sp)%x(i, j, k, 1) = rotmatrix(1, 1) * dx &
2698 + rotmatrix(1, 2) * dy &
2699 + rotmatrix(1, 3) * dz &
2701 flowdoms(mm, level, sp)%x(i, j, k, 2) = rotmatrix(2, 1) * dx &
2702 + rotmatrix(2, 2) * dy &
2703 + rotmatrix(2, 3) * dz &
2705 flowdoms(mm, level, sp)%x(i, j, k, 3) = rotmatrix(3, 1) * dx &
2706 + rotmatrix(3, 2) * dy &
2707 + rotmatrix(3, 3) * dz &
2727 integer(kind=intType),
intent(in) :: level
2731 integer :: size, procID, ierr, index
2732 integer,
dimension(mpi_status_size) :: mpiStatus
2734 integer(kind=intType) :: i, j, ii, jj, mm, idim
2735 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
2778 call mpi_isend(
recvbuffer(ii),
size, adflow_real, procid, &
2781 call echk(ierr, __file__, __line__)
2802 call mpi_irecv(
sendbuffer(ii),
size, adflow_real, procid, &
2804 call echk(ierr, __file__, __line__)
2831 flowdomsd(d1, level, mm)%x(i1, j1, k1, idim) =
flowdomsd(d1, level, mm)%x(i1, j1, k1, idim) + &
2832 flowdomsd(d2, level, mm)%x(i2, j2, k2, idim)
2853 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
2854 call echk(ierr, __file__, __line__)
2872 flowdomsd(d2, level, mm)%x(i2, j2, k2, idim) =
flowdomsd(d2, level, mm)%x(i2, j2, k2, idim) + &
2879 end do completesends
2892 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
2893 call echk(ierr, __file__, __line__)
2912 integer(kind=intType),
intent(in) :: level
2916 integer :: size, procID, ierr, index
2917 integer,
dimension(mpi_status_size) :: mpiStatus
2919 integer(kind=intType) :: i, j, ii, jj, mm
2920 integer(kind=intType) :: d1, i1, j1, k1, d2, i2, j2, k2
2964 call mpi_isend(
sendbuffer(ii),
size, adflow_real, procid, &
2967 call echk(ierr, __file__, __line__)
2988 call mpi_irecv(
recvbuffer(ii),
size, adflow_real, procid, &
2990 call echk(ierr, __file__, __line__)
3015 flowdomsd(d2, level, mm)%x(i2, j2, k2, 1) = &
3016 flowdomsd(d1, level, mm)%x(i1, j1, k1, 1)
3017 flowdomsd(d2, level, mm)%x(i2, j2, k2, 2) = &
3018 flowdomsd(d1, level, mm)%x(i1, j1, k1, 2)
3019 flowdomsd(d2, level, mm)%x(i2, j2, k2, 3) = &
3020 flowdomsd(d1, level, mm)%x(i1, j1, k1, 3)
3040 call mpi_waitany(
size,
recvrequests, index, mpistatus, ierr)
3041 call echk(ierr, __file__, __line__)
3066 end do completerecvs
3079 call mpi_waitany(
size,
sendrequests, index, mpistatus, ierr)
3080 call echk(ierr, __file__, __line__)
3105 #include <petsc/finclude/petsc.h>
3110 logical,
intent(in) :: isInflow
3111 real(kind=realtype),
dimension(:, :) :: vars
3112 integer(kind=intType) :: sps
3115 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3116 real(kind=realtype),
dimension(:),
pointer :: localptr
3134 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3135 call echk(ierr, __file__, __line__)
3138 domainloop:
do nn = 1, ndom
3140 bocoloop:
do mm = 1,
nbocos
3156 localptr(ii) =
eighth * ( &
3157 ww1(i, j, ivar) +
ww1(i + 1, j, ivar) + &
3158 ww1(i, j + 1, ivar) +
ww1(i + 1, j + 1, ivar) + &
3159 ww2(i, j, ivar) +
ww2(i + 1, j, ivar) + &
3160 ww2(i, j + 1, ivar) +
ww2(i + 1, j + 1, ivar))
3162 localptr(ii) =
eighth * ( &
3163 pp1(i, j) +
pp1(i + 1, j) + &
3164 pp1(i, j + 1) +
pp1(i + 1, j + 1) + &
3165 pp2(i, j) +
pp2(i + 1, j) + &
3166 pp2(i, j + 1) +
pp2(i + 1, j + 1))
3168 localptr(ii) =
eighth * ( &
3175 localptr(ii) =
fourth * ( &
3191 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3192 call echk(ierr, __file__, __line__)
3194 call vecplacearray(zipper%localVal, vars(:, ivar), ierr)
3195 call echk(ierr, __file__, __line__)
3198 call vecscatterbegin(zipper%scatter, exch%nodeValLocal, &
3199 zipper%localVal, insert_values, scatter_forward, ierr)
3200 call echk(ierr, __file__, __line__)
3202 call vecscatterend(zipper%scatter, exch%nodeValLocal, &
3203 zipper%localVal, insert_values, scatter_forward, ierr)
3204 call echk(ierr, __file__, __line__)
3207 call vecresetarray(zipper%localVal, ierr)
3208 call echk(ierr, __file__, __line__)
3223 #include <petsc/finclude/petsc.h>
3227 logical,
intent(in) :: isInflow
3228 real(kind=realtype),
dimension(:, :) :: vars, varsd
3229 integer(kind=intType) :: sps
3232 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3233 real(kind=realtype),
dimension(:),
pointer :: localptr
3255 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3256 call echk(ierr, __file__, __line__)
3259 domainloop:
do nn = 1, ndom
3261 bocoloop:
do mm = 1,
nbocos
3276 localptr(ii) =
eighth * ( &
3277 ww1d(i, j, ivar) +
ww1d(i + 1, j, ivar) + &
3278 ww1d(i, j + 1, ivar) +
ww1d(i + 1, j + 1, ivar) + &
3279 ww2d(i, j, ivar) +
ww2d(i + 1, j, ivar) + &
3280 ww2d(i, j + 1, ivar) +
ww2d(i + 1, j + 1, ivar))
3282 localptr(ii) =
eighth * ( &
3284 pp1d(i, j + 1) +
pp1d(i + 1, j + 1) + &
3286 pp2d(i, j + 1) +
pp2d(i + 1, j + 1))
3292 localptr(ii) =
fourth * ( &
3308 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3309 call echk(ierr, __file__, __line__)
3311 call vecplacearray(zipper%localVal, varsd(:, ivar), ierr)
3312 call echk(ierr, __file__, __line__)
3315 call vecscatterbegin(zipper%scatter, exch%nodeValLocal, &
3316 zipper%localVal, insert_values, scatter_forward, ierr)
3317 call echk(ierr, __file__, __line__)
3319 call vecscatterend(zipper%scatter, exch%nodeValLocal, &
3320 zipper%localVal, insert_values, scatter_forward, ierr)
3321 call echk(ierr, __file__, __line__)
3324 call vecresetarray(zipper%localVal, ierr)
3325 call echk(ierr, __file__, __line__)
3340 #include <petsc/finclude/petsc.h>
3345 logical,
intent(in) :: isInflow
3346 real(kind=realtype),
dimension(:, :) :: vars, varsd
3347 integer(kind=intType) :: sps
3350 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3351 real(kind=realtype),
dimension(:),
pointer :: localptr
3352 real(kind=realtype) :: tmp
3369 call vecset(exch%nodeValLocal,
zero, ierr)
3370 call echk(ierr, __file__, __line__)
3372 call vecplacearray(zipper%localVal, varsd(:, ivar), ierr)
3373 call echk(ierr, __file__, __line__)
3376 call vecscatterbegin(zipper%scatter, zipper%localVal, &
3377 exch%nodeValLocal, add_values, scatter_reverse, ierr)
3378 call echk(ierr, __file__, __line__)
3380 call vecscatterend(zipper%scatter, zipper%localVal, &
3381 exch%nodeValLocal, add_values, scatter_reverse, ierr)
3382 call echk(ierr, __file__, __line__)
3385 call vecresetarray(zipper%localVal, ierr)
3386 call echk(ierr, __file__, __line__)
3389 varsd(:, ivar) =
zero
3394 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3395 call echk(ierr, __file__, __line__)
3398 domainloop:
do nn = 1, ndom
3400 bocoloop:
do mm = 1,
nbocos
3415 tmp =
eighth * localptr(ii)
3416 ww1d(i, j, ivar) =
ww1d(i, j, ivar) + tmp
3417 ww1d(i + 1, j, ivar) =
ww1d(i + 1, j, ivar) + tmp
3418 ww1d(i, j + 1, ivar) =
ww1d(i, j + 1, ivar) + tmp
3419 ww1d(i + 1, j + 1, ivar) =
ww1d(i + 1, j + 1, ivar) + tmp
3421 ww2d(i, j, ivar) =
ww2d(i, j, ivar) + tmp
3422 ww2d(i + 1, j, ivar) =
ww2d(i + 1, j, ivar) + tmp
3423 ww2d(i, j + 1, ivar) =
ww2d(i, j + 1, ivar) + tmp
3424 ww2d(i + 1, j + 1, ivar) =
ww2d(i + 1, j + 1, ivar) + tmp
3426 tmp =
eighth * localptr(ii)
3428 pp1d(i + 1, j) =
pp1d(i + 1, j) + tmp
3429 pp1d(i, j + 1) =
pp1d(i, j + 1) + tmp
3430 pp1d(i + 1, j + 1) =
pp1d(i + 1, j + 1) + tmp
3433 pp2d(i + 1, j) =
pp2d(i + 1, j) + tmp
3434 pp2d(i, j + 1) =
pp2d(i, j + 1) + tmp
3435 pp2d(i + 1, j + 1) =
pp2d(i + 1, j + 1) + tmp
3442 tmp =
fourth * localptr(ii)
3461 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3462 call echk(ierr, __file__, __line__)
3481 #include <petsc/finclude/petsc.h>
3485 real(kind=realtype),
dimension(:, :) :: vars
3486 integer(kind=intType) :: sps
3489 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3490 real(kind=realtype),
dimension(:),
pointer :: localptr
3502 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3503 call echk(ierr, __file__, __line__)
3506 domainloop:
do nn = 1, ndom
3508 bocoloop:
do mm = 1,
nbocos
3520 localptr(ii) =
bcdata(mm)%Tp(i, j, ivar)
3539 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3540 call echk(ierr, __file__, __line__)
3542 call vecplacearray(zipper%localVal, vars(:, ivar), ierr)
3543 call echk(ierr, __file__, __line__)
3546 call vecscatterbegin(zipper%scatter, exch%nodeValLocal, &
3547 zipper%localVal, insert_values, scatter_forward, ierr)
3548 call echk(ierr, __file__, __line__)
3550 call vecscatterend(zipper%scatter, exch%nodeValLocal, &
3551 zipper%localVal, insert_values, scatter_forward, ierr)
3552 call echk(ierr, __file__, __line__)
3555 call vecresetarray(zipper%localVal, ierr)
3556 call echk(ierr, __file__, __line__)
3571 #include <petsc/finclude/petsc.h>
3576 real(kind=realtype),
dimension(:, :) :: vars, varsd
3577 integer(kind=intType) :: sps
3580 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3581 real(kind=realtype),
dimension(:),
pointer :: localptr
3596 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3597 call echk(ierr, __file__, __line__)
3600 domainloop:
do nn = 1, ndom
3602 bocoloop:
do mm = 1,
nbocos
3614 localptr(ii) =
bcdatad(mm)%Tp(i, j, ivar)
3633 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3634 call echk(ierr, __file__, __line__)
3636 call vecplacearray(zipper%localVal, varsd(:, ivar), ierr)
3637 call echk(ierr, __file__, __line__)
3640 call vecscatterbegin(zipper%scatter, exch%nodeValLocal, &
3641 zipper%localVal, insert_values, scatter_forward, ierr)
3642 call echk(ierr, __file__, __line__)
3644 call vecscatterend(zipper%scatter, exch%nodeValLocal, &
3645 zipper%localVal, insert_values, scatter_forward, ierr)
3646 call echk(ierr, __file__, __line__)
3649 call vecresetarray(zipper%localVal, ierr)
3650 call echk(ierr, __file__, __line__)
3665 #include <petsc/finclude/petsc.h>
3670 real(kind=realtype),
dimension(:, :) :: vars, varsd
3671 integer(kind=intType) :: sps
3674 integer(kind=intType) :: ii, iVar, i, j, iBeg, iEnd, jBeg, jEnd, ierr, nn, mm
3675 real(kind=realtype),
dimension(:),
pointer :: localptr
3687 call vecset(exch%nodeValLocal,
zero, ierr)
3688 call echk(ierr, __file__, __line__)
3690 call vecplacearray(zipper%localVal, varsd(:, ivar), ierr)
3691 call echk(ierr, __file__, __line__)
3694 call vecscatterbegin(zipper%scatter, zipper%localVal, &
3695 exch%nodeValLocal, add_values, scatter_reverse, ierr)
3696 call echk(ierr, __file__, __line__)
3698 call vecscatterend(zipper%scatter, zipper%localVal, &
3699 exch%nodeValLocal, add_values, scatter_reverse, ierr)
3700 call echk(ierr, __file__, __line__)
3703 call vecresetarray(zipper%localVal, ierr)
3704 call echk(ierr, __file__, __line__)
3707 varsd(:, ivar) =
zero
3712 call vecgetarrayf90(exch%nodeValLocal, localptr, ierr)
3713 call echk(ierr, __file__, __line__)
3716 domainloop:
do nn = 1, ndom
3718 bocoloop:
do mm = 1,
nbocos
3729 bcdatad(mm)%Tp(i, j, ivar) = localptr(ii)
3749 call vecrestorearrayf90(exch%nodeValLocal, localptr, ierr)
3750 call echk(ierr, __file__, __line__)
subroutine computenodaltractions(sps)
subroutine computenodaltractions_b(sps)
subroutine computenodaltractions_d(sps)
real(kind=realtype), dimension(:, :), pointer gamma2
real(kind=realtype), dimension(:, :), pointer pp1d
real(kind=realtype), dimension(:, :, :), pointer ww1d
real(kind=realtype), dimension(:, :, :), pointer ww2d
real(kind=realtype), dimension(:, :), pointer pp2d
real(kind=realtype), dimension(:, :), pointer pp1
real(kind=realtype), dimension(:, :), pointer sfaced
real(kind=realtype), dimension(:, :), pointer sface
real(kind=realtype), dimension(:, :, :), pointer ww1
real(kind=realtype), dimension(:, :), pointer pp2
real(kind=realtype), dimension(:, :, :), pointer xx
real(kind=realtype), dimension(:, :), pointer gamma1
real(kind=realtype), dimension(:, :, :), pointer ww2
real(kind=realtype), dimension(:, :, :), pointer xxd
integer(kind=inttype) ndom
type(blocktype), dimension(:, :, :), allocatable, target flowdomsd
type(blocktype), dimension(:, :, :), allocatable, target flowdoms
integer(kind=inttype), dimension(:, :), pointer orphans
real(kind=realtype), dimension(:, :, :), pointer gamma
logical addgridvelocities
integer(kind=inttype) norphans
real(kind=realtype), dimension(:, :, :), pointer p
real(kind=realtype), dimension(:, :, :, :), pointer w
type(bcdatatype), dimension(:), pointer bcdatad
integer(kind=inttype), dimension(:, :, :), pointer iblank
real(kind=realtype), dimension(:, :, :), pointer rlv
integer(kind=inttype), dimension(:), pointer bcfaceid
integer(kind=inttype) nbocos
real(kind=realtype), dimension(:, :, :), pointer rev
integer(kind=inttype), dimension(:), pointer bctype
real(kind=realtype), dimension(:), allocatable recvbuffer
integer, dimension(:), allocatable recvrequests
real(kind=realtype), dimension(:), allocatable sendbuffer
integer, dimension(:), allocatable sendrequests
type(commtype), dimension(:), allocatable commpatternnode_1st
type(internalcommtype), dimension(:), allocatable internalnode_1st
integer adflow_comm_world
integer(kind=inttype), parameter cptempcurvefits
real(kind=realtype), parameter zero
integer(kind=inttype), parameter izippflowy
integer(kind=inttype), parameter supersonicinflow
integer(kind=inttype), parameter izippwalltpy
integer(kind=inttype), parameter supersonicoutflow
real(kind=realtype), parameter eighth
integer(kind=inttype), parameter izippwallz
integer(kind=inttype), parameter izippflowx
integer(kind=inttype), parameter izippwalltvy
integer(kind=inttype), parameter subsonicoutflow
integer(kind=inttype), parameter ibcgroupoutflow
integer(kind=inttype), parameter izippflowsface
integer(kind=inttype), parameter izippwally
integer(kind=inttype), parameter izippwalltpz
integer(kind=inttype), parameter ibcgroupinflow
real(kind=realtype), parameter fourth
integer(kind=inttype), parameter izippwalltvx
integer(kind=inttype), parameter izippflowgamma
integer(kind=inttype), parameter izippwalltvz
integer(kind=inttype), parameter subsonicinflow
integer(kind=inttype), parameter nzippflowcomm
integer(kind=inttype), parameter ibcgroupwalls
integer(kind=inttype), parameter nzippwallcomm
integer(kind=inttype), parameter izippwalltpx
integer(kind=inttype), parameter izippflowp
integer(kind=inttype), parameter izippwallx
integer(kind=inttype), parameter izippflowz
subroutine computeetotblock_b(istart, iend, jstart, jend, kstart, kend, correctfork)
subroutine computeetotblock_d(istart, iend, jstart, jend, kstart, kend, correctfork)
subroutine computeetotblock(iStart, iEnd, jStart, jEnd, kStart, kEnd, correctForK)
real(kind=realtype) gammainf
real(kind=realtype) pinfcorr
real(kind=realtype), dimension(:), allocatable winf
real(kind=realtype) muinf
subroutine wallintegrationzippercomm(vars, sps)
subroutine wallintegrationzippercomm_b(vars, varsd, sps)
subroutine woverset_b(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine whalo2(level, start, end, commPressure, commGamma, commViscous)
subroutine correctperiodiccoor(level, sp, nPeriodic, periodicData)
subroutine setcommpointers(start, end, commPressure, commVarGamma, commLamVis, commEddyVis, level, sps, derivPointers, nVar, varOffset)
subroutine woverset(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine whalo2_d(level, start, end, commPressure, commGamma, commViscous)
subroutine woversetgeneric_d(nVar, level, sps, commPattern, Internal)
subroutine flowintegrationzippercomm_b(isInflow, vars, varsd, sps)
subroutine whalo1to1intgeneric_b(nVar, level, sps, commPattern, internal)
subroutine whalo1to1realgeneric(nVar, level, sps, commPattern, internal)
subroutine orphanaverage(wstart, wend, calcPressure, calcGamma, calcLamVis, calcEddyVis)
subroutine whalo1to1intgeneric(nVar, level, sps, commPattern, internal)
subroutine woverset_d(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine whalo1to1_d(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine wallintegrationzippercomm_d(vars, varsd, sps)
subroutine whalo1to1_b(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine whalo1(level, start, end, commPressure, commGamma, commViscous)
subroutine whalo1to1realgeneric_b(nVar, level, sps, commPattern, internal)
subroutine woversetgeneric(nVar, level, sps, commPattern, Internal)
subroutine flowintegrationzippercomm_d(isInflow, vars, varsd, sps)
subroutine correctperiodicvelocity(level, sp, nPeriodic, periodicData)
subroutine woversetgeneric_b(nVar, level, sps, commPattern, Internal)
subroutine exchangecoor_d(level)
subroutine flowintegrationzippercomm(isInflow, vars, sps)
subroutine whalo2_b(level, start, end, commPressure, commGamma, commViscous)
subroutine exchangecoor(level)
subroutine reshalo1(level, start, end)
subroutine whalo1to1(level, start, end, commPressure, commVarGamma, commLamVis, commEddyVis, commPattern, internal)
subroutine exchangecoor_b(level)
type(zippermesh), dimension(nfamexchange), target zippermeshes
type(familyexchange), dimension(:, :), allocatable, target bcfamexchange
logical function iswalltype(bType)
subroutine setbcpointers(nn, spatialPointers)
subroutine setpointers_d(nn, level, sps)
subroutine setbcpointers_d(nn, spatialpointers)
subroutine echk(errorcode, file, line)
subroutine setpointers_b(nn, level, sps)
logical function getcorrectfork()
subroutine setpointers(nn, mm, ll)