65 CHARACTER(CHRLEN) :: RCSIdentString = &
66 '$RCSfile: RFLO_ModMoveGridElliptFra.F90,v $ $Revision: 1.16 $'
110 include
'roccomf90.h'
114 TYPE(t_region
),
POINTER :: regions(:)
117 INTEGER :: ireg, iter, ipatch, ijk
124 REAL(RFREAL) :: resid, globalresid
125 REAL(RFREAL),
POINTER :: xyz(:,:), xyzold(:,:)
131 DOUBLE PRECISION :: dalpha
136 global => regions(1)%global
139 'RFLO_ModMoveGridElliptFra.F90' )
144 dalpha = global%dtMin/global%dTimeSystem
145 CALL com_call_function( global%genxHandleGm,1,dalpha )
160 DO ireg=1,global%nRegions
161 IF (regions(ireg)%procid==global%myProcid .AND. &
162 regions(ireg)%active==active .AND. &
163 regions(ireg)%mixtInput%moveGrid)
THEN
167 grid => regions(ireg)%levels(1)%grid
168 gridold => regions(ireg)%levels(1)%gridOld
172 gridold%xyzOld,
grid%xyz )
185 IF (global%moveGridNiter < 1 .AND. &
186 global%moveGridViter < 1 .AND. &
187 global%moveGridSiter < 1)
THEN
188 IF (global%verbLevel >= verbose_high)
THEN
189 IF (global%myProcid == masterproc)
THEN
190 WRITE(stdout,4000) solver_name,global%skewness,global%minVol
191 WRITE(stdout,1000) solver_name, &
192 global%moveGridNiter, global%moveGridAmplifX, &
193 global%moveGridAmplifY, global%moveGridAmplifZ, &
194 global%moveGridPower, global%moveGridOrthDir, &
195 global%moveGridOrthWghtX, global%moveGridOrthWghtY, &
196 global%moveGridOrthWghtZ
204 DO iter=1,global%moveGridNiter
211 IF (global%verbLevel >= verbose_high)
THEN
213 CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
214 masterproc,global%mpiComm,global%mpierr )
215 IF (global%mpierr /= 0) CALL
errorstop( global,err_mpi_trouble,&
220 IF (global%myProcid == masterproc)
THEN
221 WRITE(stdout,4000) solver_name,global%skewness,global%minVol
222 WRITE(stdout,2000) solver_name, &
223 global%moveGridNiter, global%moveGridAmplifX, &
224 global%moveGridAmplifY,global%moveGridAmplifZ, &
225 global%moveGridPower,global%moveGridOrthDir, &
226 global%moveGridOrthWghtX, global%moveGridOrthWghtY, &
227 global%moveGridOrthWghtZ,
sqrt(globalresid)
234 DO ireg=1,global%nRegions
235 IF (regions(ireg)%procid==global%myProcid .AND. &
236 regions(ireg)%active==active .AND. &
237 regions(ireg)%mixtInput%moveGrid)
THEN
241 xyz => regions(ireg)%levels(1)%grid%xyz
242 xyzold => regions(ireg)%levels(1)%gridOld%xyz
244 DO ijk=lbound(xyz,2),ubound(xyz,2)
245 xyz(xcoord,ijk) = xyz(xcoord,ijk) - xyzold(xcoord,ijk)
246 xyz(ycoord,ijk) = xyz(ycoord,ijk) - xyzold(ycoord,ijk)
247 xyz(zcoord,ijk) = xyz(zcoord,ijk) - xyzold(zcoord,ijk)
252 grid => regions(ireg)%levels(1)%grid
253 gridold => regions(ireg)%levels(1)%gridOld
254 grid%boundMoved(:) = .true.
255 grid%edgeMoved(:) = .true.
257 DO ipatch=1,regions(ireg)%nPatches
258 patch => regions(ireg)%levels(1)%patches(ipatch)
259 bctype =
patch%bcType
261 IF (bctype .eq. bc_symmetry)
THEN
271 gridold%xyzOld,
grid%xyz )
280 gridold%xyzOld,
grid%xyz )
284 IF (regions(ireg)%blockShape==region_shape_normal)
THEN
285 DO iter=1,global%moveGridViter
302 IF (global%verbLevel >= verbose_high)
THEN
304 CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
305 masterproc,global%mpiComm,global%mpierr )
306 IF (global%mpierr /= 0) CALL
errorstop( global,err_mpi_trouble,&
311 IF (global%myProcid == masterproc)
THEN
312 WRITE(stdout,3000) solver_name,global%moveGridViter,
sqrt(globalresid)
320 DO ireg=1,global%nRegions
321 IF (regions(ireg)%procid==global%myProcid .AND. &
322 regions(ireg)%active==active .AND. &
323 regions(ireg)%mixtInput%moveGrid)
THEN
327 IF (regions(ireg)%mixtInput%faceEdgeAvg==fe_avg_linear) &
343 1000
FORMAT(
a,1
x,
'Global-Elliptic-PDE grid motion:', &
344 i5,1
x,4(1pe9.2),i4,3(1pe9.2))
345 2000
FORMAT(
a,1
x,
'Global-Elliptic-PDE grid motion:', &
346 i5,1
x,4(1pe9.2),i4,3(1pe9.2),1pe12.4)
347 3000
FORMAT(
a,1
x,
'Global-Elliptic-PDE grid motion:', &
349 4000
FORMAT(
a,1
x,
'global skewness, minvol:',2(1pe14.5))
385 TYPE(t_region
),
POINTER :: regions(:)
388 INTEGER :: ireg, iter
396 global => regions(1)%global
399 'RFLO_ModMoveGridElliptFra.F90' )
405 DO ireg=1,global%nRegions
406 IF (regions(ireg)%procid==global%myProcid .AND. &
407 regions(ireg)%active==active .AND. &
408 regions(ireg)%mixtInput%moveGrid)
THEN
410 grid => regions(ireg)%levels(1)%grid
411 gridold => regions(ireg)%levels(1)%gridOld
416 gridold%indSvel =
grid%indSvel
417 gridold%ipc =
grid%ipc
418 gridold%jpc =
grid%jpc
419 gridold%kpc =
grid%kpc
420 gridold%xyz(:,:) =
grid%xyz(:,:)
421 gridold%si(:,:) =
grid%si(:,:)
422 gridold%sj(:,:) =
grid%sj(:,:)
423 gridold%sk(:,:) =
grid%sk(:,:)
424 gridold%vol(:) =
grid%vol(:)
456 DO ireg=1,global%nRegions
457 IF (regions(ireg)%procid==global%myProcid .AND. &
458 regions(ireg)%active==active .AND. &
459 regions(ireg)%mixtInput%moveGrid)
THEN
461 grid => regions(ireg)%levels(1)%grid
462 gridold => regions(ireg)%levels(1)%gridOld
466 grid%edgeMoved(:) = .false.
469 grid%arcLen34,
grid%arcLen56,gridold%xyzOld, &
470 gridold%xyz,
grid%xyz )
474 IF (global%moveGridNiter < 1)
THEN
478 gridold%xyzOld,gridold%xyz,
grid%xyz )
486 gridold%xyzOld,
grid%xyz )
522 arclen12,arclen34,arclen56,xyzref,xyzold, &
529 #include "Indexing.h"
532 LOGICAL :: boundmoved(6), allexternal(6), edgemoved(12)
534 REAL(RFREAL),
POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
535 REAL(RFREAL),
POINTER :: dnode(:,:), xyzold(:,:), xyzref(:,:)
537 TYPE(t_region
) :: region
540 INTEGER :: iedge, ind
544 INTEGER :: indbeg, indend, ijkn, ijkn1, ijknbeg, ijknend, inoff, ijnoff
545 INTEGER :: switch(12,11), intertype, iedgeglo
547 REAL(RFREAL) :: arclen, ds,
s, dn(3), dnbeg(3), dnend(3)
553 'RFLO_ModMoveGridElliptFra.F90')
575 switch( 1,:) = (/5, 6, 1, 3, 56,
kpnbeg, kpnend,
ipnbeg,
jpnbeg, 1, 2/)
576 switch( 2,:) = (/3, 4, 1, 6, 34,
jpnbeg,
jpnend, kpnend,
ipnbeg, 2, 3/)
577 switch( 3,:) = (/5, 6, 1, 4, 56,
kpnbeg, kpnend,
ipnbeg,
jpnend, 4, 3/)
578 switch( 4,:) = (/3, 4, 1, 5, 34,
jpnbeg,
jpnend,
kpnbeg,
ipnbeg, 1, 4/)
579 switch( 5,:) = (/5, 6, 2, 3, 56,
kpnbeg, kpnend,
ipnend,
jpnbeg, 5, 6/)
580 switch( 6,:) = (/3, 4, 2, 6, 34,
jpnbeg,
jpnend, kpnend,
ipnend, 6, 7/)
581 switch( 7,:) = (/5, 6, 2, 4, 56,
kpnbeg, kpnend,
ipnend,
jpnend, 8, 7/)
582 switch( 8,:) = (/3, 4, 2, 5, 34,
jpnbeg,
jpnend,
kpnbeg,
ipnend, 5, 8/)
583 switch( 9,:) = (/1, 2, 3, 5, 12,
ipnbeg,
ipnend,
jpnbeg,
kpnbeg, 1, 5/)
584 switch(10,:) = (/1, 2, 3, 6, 12,
ipnbeg,
ipnend,
jpnbeg, kpnend, 2, 6/)
585 switch(11,:) = (/1, 2, 4, 5, 12,
ipnbeg,
ipnend,
jpnend,
kpnbeg, 4, 8/)
586 switch(12,:) = (/1, 2, 4, 6, 12,
ipnbeg,
ipnend,
jpnend, kpnend, 3, 7/)
590 IF (boundmoved(1))
THEN
591 edgemoved( 1) = .true.; edgemoved( 2) = .true.
592 edgemoved( 3) = .true.; edgemoved( 4) = .true.
594 IF (boundmoved(2))
THEN
595 edgemoved( 5) = .true.; edgemoved( 6) = .true.
596 edgemoved( 7) = .true.; edgemoved( 8) = .true.
598 IF (boundmoved(3))
THEN
599 edgemoved( 1) = .true.; edgemoved( 5) = .true.
600 edgemoved( 9) = .true.; edgemoved(10) = .true.
602 IF (boundmoved(4))
THEN
603 edgemoved( 3) = .true.; edgemoved( 7) = .true.
604 edgemoved(11) = .true.; edgemoved(12) = .true.
606 IF (boundmoved(5))
THEN
607 edgemoved( 4) = .true.; edgemoved( 8) = .true.
608 edgemoved( 9) = .true.; edgemoved(11) = .true.
610 IF (boundmoved(6))
THEN
611 edgemoved( 2) = .true.; edgemoved( 6) = .true.
612 edgemoved(10) = .true.; edgemoved(12) = .true.
618 IF (.NOT.edgemoved(iedge))
THEN
620 edgemoved(iedge) = .true.
623 indbeg = switch(iedge,6)
624 indend = switch(iedge,7)
625 l1c = switch(iedge,8)
626 l2c = switch(iedge,9)
629 IF (iedge==11) iedgeglo=12
630 IF (iedge==12) iedgeglo=11
631 interact = region%levels(ilev)%edgeCells(iedgeglo)%interact
632 intertype = region%levels(ilev)%edgeCells(iedgeglo)%interType
634 IF (((region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==1 .OR. &
635 region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==1) .AND. &
636 ((interact .EQV. .true.) .AND. (intertype==edge_interact_full))) &
638 region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==2 .OR. &
639 region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==2)
THEN
646 DO ind=indbeg+1,indend-1
647 IF (switch(iedge,5) == 12)
THEN
648 ijkn = indijk(ind ,l1c,l2c,inoff,ijnoff)
649 ijkn1 = indijk(ind-1 ,l1c,l2c,inoff,ijnoff)
650 ijknbeg = indijk(indbeg,l1c,l2c,inoff,ijnoff)
651 ijknend = indijk(indend,l1c,l2c,inoff,ijnoff)
652 arclen = arclen12(l1c,l2c)
653 ELSE IF (switch(iedge,5) == 34)
THEN
654 ijkn = indijk(l2c,ind ,l1c,inoff,ijnoff)
655 ijkn1 = indijk(l2c,ind-1 ,l1c,inoff,ijnoff)
656 ijknbeg = indijk(l2c,indbeg,l1c,inoff,ijnoff)
657 ijknend = indijk(l2c,indend,l1c,inoff,ijnoff)
658 arclen = arclen34(l1c,l2c)
659 ELSE IF (switch(iedge,5) == 56)
THEN
660 ijkn = indijk(l1c,l2c,ind ,inoff,ijnoff)
661 ijkn1 = indijk(l1c,l2c,ind-1 ,inoff,ijnoff)
662 ijknbeg = indijk(l1c,l2c,indbeg,inoff,ijnoff)
663 ijknend = indijk(l1c,l2c,indend,inoff,ijnoff)
664 arclen = arclen56(l1c,l2c)
666 dnbeg(:) = dnode(:,ijknbeg) + xyzold(:,ijknbeg) - xyzref(:,ijknbeg)
667 dnend(:) = dnode(:,ijknend) + xyzold(:,ijknend) - xyzref(:,ijknend)
669 ds = ds +
sqrt((xyzref(xcoord,ijkn)-xyzref(xcoord,ijkn1))**2 + &
670 (xyzref(ycoord,ijkn)-xyzref(ycoord,ijkn1))**2 + &
671 (xyzref(zcoord,ijkn)-xyzref(zcoord,ijkn1))**2)
675 dnode(:,ijkn) = dn(:) + xyzref(:,ijkn) - xyzold(:,ijkn)
710 arclen12,arclen34,arclen56, &
717 #include "Indexing.h"
720 LOGICAL :: boundmoved(6), edgemoved(12)
722 REAL(RFREAL),
POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
723 REAL(RFREAL),
POINTER :: dnode(:,:), xyzold(:,:)
725 TYPE(t_region
) :: region
728 INTEGER :: ibound, l1, l2
732 INTEGER :: l1b, l1e, l2b, l2e, lc, ijkn, ijke(4), ijkem(4), inoff, ijnoff
733 INTEGER :: switch(6,9)
737 REAL(RFREAL) :: arclen(4), ds(4),
s(4)
738 REAL(RFREAL) :: corner(3,8), e1(3), e2(3), e3(3), e4(3), &
739 p1(3), p2(3), p3(3), p4(3), dn(3)
744 'RFLO_ModMoveGridElliptFra.F90' )
769 corner(:,2) = dnode(:,indijk(
ipnbeg,
jpnbeg,kpnend,inoff,ijnoff))
770 corner(:,3) = dnode(:,indijk(
ipnbeg,
jpnend,kpnend,inoff,ijnoff))
773 corner(:,6) = dnode(:,indijk(
ipnend,
jpnbeg,kpnend,inoff,ijnoff))
774 corner(:,7) = dnode(:,indijk(
ipnend,
jpnend,kpnend,inoff,ijnoff))
780 IF ((.NOT.boundmoved(ibound)) .AND. &
781 (edgemoved(switch(ibound,1)) .OR. edgemoved(switch(ibound,2)) .OR. &
782 edgemoved(switch(ibound,3)) .OR. edgemoved(switch(ibound,4))))
THEN
787 l1b = switch(ibound,5)
788 l1e = switch(ibound,6)
789 l2b = switch(ibound,7)
790 l2e = switch(ibound,8)
791 lc = switch(ibound,9)
793 IF (ibound == 1)
THEN
798 ELSE IF (ibound == 2)
THEN
803 ELSE IF (ibound == 3)
THEN
808 ELSE IF (ibound == 4)
THEN
813 ELSE IF (ibound == 5)
THEN
818 ELSE IF (ibound == 6)
THEN
831 IF (ibound==1 .OR. ibound==2)
THEN
832 ijkn = indijk(lc,l1 ,l2 ,inoff,ijnoff)
833 ijke(1) = indijk(lc,
jpnbeg,l2 ,inoff,ijnoff)
834 ijkem(1) = indijk(lc,
jpnbeg,l2-1 ,inoff,ijnoff)
835 ijke(2) = indijk(lc,
jpnend,l2 ,inoff,ijnoff)
836 ijkem(2) = indijk(lc,
jpnend,l2-1 ,inoff,ijnoff)
837 ijke(3) = indijk(lc,l1 ,
kpnbeg,inoff,ijnoff)
838 ijkem(3) = indijk(lc,l1-1 ,
kpnbeg,inoff,ijnoff)
839 ijke(4) = indijk(lc,l1 ,kpnend,inoff,ijnoff)
840 ijkem(4) = indijk(lc,l1-1 ,kpnend,inoff,ijnoff)
841 arclen(1) = arclen56(lc,
jpnbeg)
842 arclen(2) = arclen56(lc,
jpnend)
843 arclen(3) = arclen34(
kpnbeg,lc)
844 arclen(4) = arclen34(kpnend,lc)
845 ELSE IF (ibound==3 .OR. ibound==4)
THEN
846 ijkn = indijk(l2 ,lc,l1 ,inoff,ijnoff)
847 ijke(1) = indijk(l2 ,lc,
kpnbeg,inoff,ijnoff)
848 ijkem(1) = indijk(l2-1 ,lc,
kpnbeg,inoff,ijnoff)
849 ijke(2) = indijk(l2 ,lc,kpnend,inoff,ijnoff)
850 ijkem(2) = indijk(l2-1 ,lc,kpnend,inoff,ijnoff)
851 ijke(3) = indijk(
ipnbeg,lc,l1 ,inoff,ijnoff)
852 ijkem(3) = indijk(
ipnbeg,lc,l1-1 ,inoff,ijnoff)
853 ijke(4) = indijk(
ipnend,lc,l1 ,inoff,ijnoff)
854 ijkem(4) = indijk(
ipnend,lc,l1-1 ,inoff,ijnoff)
855 arclen(1) = arclen12(lc,
kpnbeg)
856 arclen(2) = arclen12(lc,kpnend)
857 arclen(3) = arclen56(
ipnbeg,lc)
858 arclen(4) = arclen56(
ipnend,lc)
859 ELSE IF (ibound==5 .OR. ibound==6)
THEN
860 ijkn = indijk(l1 ,l2 ,lc,inoff,ijnoff)
861 ijke(1) = indijk(
ipnbeg,l2 ,lc,inoff,ijnoff)
862 ijkem(1) = indijk(
ipnbeg,l2-1 ,lc,inoff,ijnoff)
863 ijke(2) = indijk(
ipnend,l2 ,lc,inoff,ijnoff)
864 ijkem(2) = indijk(
ipnend,l2-1 ,lc,inoff,ijnoff)
865 ijke(3) = indijk(l1 ,
jpnbeg,lc,inoff,ijnoff)
866 ijkem(3) = indijk(l1-1 ,
jpnbeg,lc,inoff,ijnoff)
867 ijke(4) = indijk(l1 ,
jpnend,lc,inoff,ijnoff)
868 ijkem(4) = indijk(l1-1 ,
jpnend,lc,inoff,ijnoff)
869 arclen(1) = arclen34(lc,
ipnbeg)
870 arclen(2) = arclen34(lc,
ipnend)
871 arclen(3) = arclen12(
jpnbeg,lc)
872 arclen(4) = arclen12(
jpnend,lc)
876 sqrt((xyzold(xcoord,ijke(1))-xyzold(xcoord,ijkem(1)))**2 + &
877 (xyzold(ycoord,ijke(1))-xyzold(ycoord,ijkem(1)))**2 + &
878 (xyzold(zcoord,ijke(1))-xyzold(zcoord,ijkem(1)))**2)
880 sqrt((xyzold(xcoord,ijke(2))-xyzold(xcoord,ijkem(2)))**2 + &
881 (xyzold(ycoord,ijke(2))-xyzold(ycoord,ijkem(2)))**2 + &
882 (xyzold(zcoord,ijke(2))-xyzold(zcoord,ijkem(2)))**2)
886 sqrt((xyzold(xcoord,ijke(3))-xyzold(xcoord,ijkem(3)))**2 + &
887 (xyzold(ycoord,ijke(3))-xyzold(ycoord,ijkem(3)))**2 + &
888 (xyzold(zcoord,ijke(3))-xyzold(zcoord,ijkem(3)))**2)
890 sqrt((xyzold(xcoord,ijke(4))-xyzold(xcoord,ijkem(4)))**2 + &
891 (xyzold(ycoord,ijke(4))-xyzold(ycoord,ijkem(4)))**2 + &
892 (xyzold(zcoord,ijke(4))-xyzold(zcoord,ijkem(4)))**2)
893 s(:) = ds(:)/arclen(:)
894 e1(:) = dnode(:,ijke(1))
895 e2(:) = dnode(:,ijke(2))
896 e3(:) = dnode(:,ijke(3))
897 e4(:) = dnode(:,ijke(4))
898 CALL
rflo_tfint2d(
s(1),
s(2),
s(3),
s(4),e1,e2,e3,e4,
p1,p2,p3,p4,dn )
899 dnode(:,ijkn) = dn(:)
937 TYPE(t_region
),
POINTER :: regions(:)
940 INTEGER :: ireg, ipatch, ipass
943 INTEGER :: bctype, iregsrc, ipatchsrc
945 TYPE(t_grid),
POINTER ::
grid, gridold, gridsrc
951 global => regions(1)%global
954 'RFLO_ModMoveGridElliptFra.F90' )
962 DO ireg=1,global%nRegions
963 IF (regions(ireg)%procid==global%myProcid .AND. &
964 regions(ireg)%active==active .AND. &
965 regions(ireg)%mixtInput%moveGrid)
THEN
967 grid => regions(ireg)%levels(1)%grid
968 gridold => regions(ireg)%levels(1)%gridOld
970 DO ipatch=1,regions(ireg)%nPatches
971 patch => regions(ireg)%levels(1)%patches(ipatch)
972 bctype =
patch%bcType
973 IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
974 (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
975 (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range))
THEN
976 iregsrc =
patch%srcRegion
977 ipatchsrc =
patch%srcPatch
978 patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
979 gridsrc => regions(iregsrc)%levels(1)%grid
981 IF (regions(iregsrc)%procid == global%myProcid)
THEN
983 patch,patchsrc,.false., &
984 grid%xyz,gridsrc%xyz )
988 grid%arcLen56,gridold%xyzOld, &
989 gridold%xyz,
grid%xyz )
993 gridold%xyzOld,
grid%xyz )
1006 DO ireg=1,global%nRegions
1007 IF (regions(ireg)%procid==global%myProcid .AND. &
1008 regions(ireg)%active==active .AND. &
1009 regions(ireg)%mixtInput%moveGrid)
THEN
1011 grid => regions(ireg)%levels(1)%grid
1012 gridold => regions(ireg)%levels(1)%gridOld
1014 DO ipatch=1,regions(ireg)%nPatches
1015 patch => regions(ireg)%levels(1)%patches(ipatch)
1016 bctype =
patch%bcType
1017 IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
1018 (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
1019 (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range))
THEN
1020 iregsrc =
patch%srcRegion
1021 ipatchsrc =
patch%srcPatch
1022 patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
1023 gridsrc => regions(iregsrc)%levels(1)%grid
1025 IF (regions(iregsrc)%procid /= global%myProcid)
THEN
1031 grid%arcLen56,gridold%xyzOld, &
1032 gridold%xyz,
grid%xyz )
1036 gridold%xyzOld,
grid%xyz )
1046 DO ireg=1,global%nRegions
1047 IF (regions(ireg)%procid==global%myProcid .AND. &
1048 regions(ireg)%active==active .AND. &
1049 regions(ireg)%mixtInput%moveGrid)
THEN
subroutine, public rflo_elliptgridsmoo(regions, resid)
subroutine rflo_copygeometrydummy(region)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jpnbeg
subroutine, public rflo_elliptgridsmooregion(region, resid)
subroutine rflo_calccellcentroids(region)
subroutine, public rflo_movegridelliptfra(regions)
subroutine rflo_arclengthbounds(region, xyz, arcLen12, arcLen34, arcLen56)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE kpnbeg
subroutine rflo_tfint2d(s1, s2, s3, s4, e1, e2, e3, e4, p1, p2, p3, p4, xyz)
subroutine rflo_c2eavgcoeffs(region)
subroutine, public rflo_mgframebroadcast(regions, iselect, iter)
subroutine registerfunction(global, funName, fileName)
subroutine rflo_mgelliptboundaries(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine rflo_exchangednoderecv(region, regionSrc, patch, patchSrc, average, dNode)
subroutine rflo_c2favgcoeffs(region)
subroutine rflo_changeinteriorgrid(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, xyz)
subroutine, public rflo_gridqualityglobal(regions)
subroutine, public rflo_mgelliptbnddeformation(region)
subroutine rflo_movegridinterfaces(regions)
subroutine rflo_calccontrolvolumes(region)
subroutine rflo_tfint1d(s, p1, p2, xyz)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jpnend
subroutine rflo_getnodeoffset(region, iLev, iNodeOffset, ijNodeOffset)
subroutine rflo_mgframeorthoshift(regions)
subroutine rflo_mgframemovecorners(regions)
subroutine rflo_calcfacevectors(region)
subroutine rflo_exchangegeometry(regions)
subroutine rflo_generatecoarsegrids(region)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ipnbeg
subroutine rflo_mgelliptsurfacesfra(regions, someMoved)
subroutine rflo_clearsendrequests(regions, iReg, geometry)
subroutine rflo_calcgridspeeds(region)
subroutine rflo_mgelliptinterfacesfra(regions)
subroutine rflo_getdimensphysnodes(region, iLev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend)
subroutine rflo_exchangednodecopy(region, regionSrc, patch, patchSrc, average, dNode, dNodeSrc)
subroutine errorstop(global, errorCode, errorLine, addMessage)
subroutine rflo_mgelliptedgesfra(region, boundMoved, allExternal, edgeMoved, arcLen12, arcLen34, arcLen56, xyzRef, xyzOld, dNode)
subroutine deregisterfunction(global)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ipnend
subroutine rflo_exchangednodesend(region, regionSrc, patch, dNode)
subroutine rflo_checkmetrics(iReg, region)
subroutine rflo_laplacegridsmoo(regions, resid)
subroutine rflo_calcfacecentroids(region)
subroutine rflo_mgframecorrectneighbors(regions)