Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLO_ModMoveGridNconform.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 ! ******************************************************************************
24 !
25 ! Purpose: Suite for frame grid-motion routines.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLO_ModMoveGridNconform.F90,v 1.37 2009/08/27 14:04:50 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2004 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY : t_global
42  USE moddatastruct, ONLY: t_region
43  USE modgrid, ONLY : t_grid
44  USE modbndpatch, ONLY : t_patch
45  USE modparameters
46  USE moddatatypes
47  USE moderror
48  USE modmpi
49 
50  IMPLICIT NONE
51 
52  PRIVATE
53  PUBLIC :: rflo_movegridframe, &
57 
63 
64 ! private : RFLO_MgFrameSurface
65 ! RFLO_MgFrameInterfaces
66 
67 ! ******************************************************************************
68 ! Declarations and definitions
69 ! ******************************************************************************
70 
71  CHARACTER(CHRLEN) :: rcsidentstring = &
72  '$RCSfile: RFLO_ModMoveGridNconform.F90,v $ $Revision: 1.37 $'
73 
74 ! ******************************************************************************
75 ! Routines
76 ! ******************************************************************************
77 
78  CONTAINS
79 
80 !******************************************************************************
81 !
82 ! Purpose: redistribute grid nodes according to the movement of the
83 ! boundaries. This function smoothes the grid globally by
84 ! volume mesh smoothing based on Laplacian propagation.
85 !
86 ! Description: none.
87 !
88 ! Input: regions = data of all grid regions.
89 !
90 ! Output: regions%levels%grid%xyz = new grid coordinates.
91 !
92 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
93 ! is applied to the finest grid first.
94 !
95 !******************************************************************************
96 
97 SUBROUTINE rflo_movegridframe( regions )
98 
107 
108  IMPLICIT NONE
109 
110 #ifdef GENX
111  include 'roccomf90.h'
112 #endif
113 
114 ! ... parameters
115  TYPE(t_region), POINTER :: regions(:)
116 
117 ! ... loop variables
118  INTEGER :: ireg, iter, ipatch, ijk
119 
120 ! ... local variables
121  LOGICAL :: somemoved, someremesh
122 
123  INTEGER :: bctype, iremesh, jremesh, nremesh, itype
124 
125  REAL(RFREAL) :: resid, globalresid
126  REAL(RFREAL), POINTER :: xyz(:,:), xyzold(:,:)
127 
128  TYPE(t_grid), POINTER :: grid, gridold
129  TYPE(t_global), POINTER :: global
130  TYPE(t_patch), POINTER :: patch
131 #ifdef GENX
132  DOUBLE PRECISION :: dalpha
133 #endif
134 
135 !******************************************************************************
136 
137  global => regions(1)%global
138 
139  CALL registerfunction( global,'RFLO_MoveGridFrame',&
140  'RFLO_ModMoveGridNconform.F90' )
141 
142  itype=2
143 
144 #ifdef GENX
145 ! update geometry buffers -----------------------------------------------------
146 
147  dalpha = global%dtMin/global%dTimeSystem
148  CALL com_call_function( global%genxHandleGm,1,dalpha )
149 #endif
150 
151 ! receive and distribute deformations for each region -------------------------
152 
153  CALL rflo_mgframesurfaces( regions,somemoved,itype )
154 
155 ! fix interfaces between regions ----------------------------------------------
156 
157  IF (somemoved) THEN
158  CALL rflo_mgframeinterfaces( regions,itype )
159  ENDIF
160 
161 ! update grid, dummy, corner and edge cells -----------------------------------
162 
163  DO ireg=1,global%nRegions
164  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
165  regions(ireg)%active==active .AND. & ! on my processor
166  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
167 
168 ! --- change the interior grid
169 
170  grid => regions(ireg)%levels(1)%grid
171  gridold => regions(ireg)%levels(1)%gridOld
172  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
173  grid%edgeMoved,grid%arcLen12, &
174  grid%arcLen34,grid%arcLen56, &
175  gridold%xyzOld,grid%xyz )
176 
177 ! --- update coarse grids and dummy cells
178 
179  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
180  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
181  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
182  ENDIF ! region on this processor and active, grid moving
183  ENDDO ! iReg
184  CALL rflo_exchangegeometry( regions ) ! exchange geometry
185 
186 ! smooth grid by solving Laplace equation -------------------------------------
187 
188  IF (global%moveGridNiter < 1) THEN
189  IF (global%verbLevel >= verbose_high) THEN
190  IF (global%myProcid == masterproc) THEN
191  WRITE(stdout,4000) solver_name,global%skewness,global%minVol
192  WRITE(stdout,1000) solver_name, global%moveGridNiter, &
193  global%moveGridAmplifX,global%moveGridAmplifY, &
194  global%moveGridAmplifZ,global%moveGridPower
195  ENDIF ! masterproc
196  ENDIF ! verbLevel
197  goto 888
198  ENDIF ! niter<1
199 
200  IF (somemoved) THEN
201  DO iter=1,global%moveGridNiter
202  CALL rflo_laplacegridsmoo( regions,resid )
203  ENDDO
204 
205  IF (global%verbLevel >= verbose_high) THEN
206 #ifdef MPI
207  CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
208  masterproc,global%mpiComm,global%mpierr )
209  IF (global%mpierr /= 0) CALL errorstop( global,err_mpi_trouble,__line__ )
210 #else
211  globalresid = resid
212 #endif
213  IF (global%myProcid == masterproc) THEN
214  WRITE(stdout,4000) solver_name,global%skewness,global%minVol
215 
216  IF (global%moveGridScheme==movegrid_frame) THEN
217  WRITE(stdout,2000) solver_name, &
218  global%moveGridNiter, &
219  global%moveGridAmplifX,global%moveGridAmplifY, &
220  global%moveGridAmplifZ,global%moveGridPower, &
221  sqrt(globalresid)
222  ELSEIF (global%moveGridScheme==movegrid_foms) THEN
223  WRITE(stdout,3000) solver_name, &
224  global%moveGridNiter, &
225  global%moveGridAmplifX,global%moveGridAmplifY, &
226  global%moveGridAmplifZ,global%moveGridPower, &
227  global%moveGridWeight,global%moveGridOrthCell, &
228  sqrt(globalresid)
229  ENDIF
230  ENDIF
231  ENDIF ! verbLevel
232  ENDIF ! someMoved
233 
234 ! update grid, dummy, corner and edge cells -----------------------------------
235 
236  DO ireg=1,global%nRegions
237  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
238  regions(ireg)%active==active .AND. & ! on my processor
239  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
240 
241 ! --- change xyz from coordinates to deformations
242 
243  xyz => regions(ireg)%levels(1)%grid%xyz
244  xyzold => regions(ireg)%levels(1)%gridOld%xyz
245 
246  DO ijk=lbound(xyz,2),ubound(xyz,2)
247  xyz(xcoord,ijk) = xyz(xcoord,ijk) - xyzold(xcoord,ijk)
248  xyz(ycoord,ijk) = xyz(ycoord,ijk) - xyzold(ycoord,ijk)
249  xyz(zcoord,ijk) = xyz(zcoord,ijk) - xyzold(zcoord,ijk)
250  ENDDO
251 
252 ! --- redistribute deformations at boundaries
253 
254  grid => regions(ireg)%levels(1)%grid
255  gridold => regions(ireg)%levels(1)%gridOld
256  grid%boundMoved(:) = .true.
257  grid%edgeMoved(:) = .true.
258  DO ipatch=1,regions(ireg)%nPatches
259  patch => regions(ireg)%levels(1)%patches(ipatch)
260  bctype = patch%bcType
261 ! IF ((bcType>=BC_SYMMETRY .AND. bcType<=BC_SYMMETRY+BC_RANGE)) THEN
262 ! grid%boundMoved(patch%lbound) = .false.
263 ! ENDIF ! bcType
264  IF ((bctype==bc_symmetry)) THEN
265  grid%boundMoved(patch%lbound) = .false.
266  ENDIF ! bcType
267  ENDDO ! iPatch
268  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
269  grid%edgeMoved,grid%arcLen12, &
270  grid%arcLen34,grid%arcLen56, &
271  gridold%xyzOld,grid%xyz )
272 
273 ! --- change xyz from deformations to coordinates
274 
275  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
276  grid%edgeMoved,grid%arcLen12, &
277  grid%arcLen34,grid%arcLen56, &
278  gridold%xyzOld,grid%xyz )
279 
280 ! --- update coarse grids and dummy cells
281 
282  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
283  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
284  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
285  ENDIF ! region on this processor and active, grid moving
286  ENDDO ! iReg
287 
288  CALL rflo_exchangegeometry( regions ) ! exchange geometry
289 
290 888 CONTINUE
291 
292 ! calculate new metrics and grid speeds ---------------------------------------
293 
294  someremesh = .false.
295  iremesh = 0
296  DO ireg=1,global%nRegions
297  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
298  regions(ireg)%active==active .AND. & ! on my processor
299  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
300  CALL rflo_calcfacevectors( regions(ireg) ) ! faces
301  CALL rflo_calccontrolvolumes( regions(ireg) ) ! volumes
302  CALL rflo_calccellcentroids( regions(ireg) ) ! cell centroids
303  IF (global%moveGridScheme==movegrid_foms) &
304  CALL rflo_calcfacecentroids( regions(ireg) ) ! face centroids
305  IF (regions(ireg)%mixtInput%faceEdgeAvg==fe_avg_linear) &
306  CALL rflo_c2favgcoeffs( regions(ireg) ) ! cell2face averaging
307  CALL rflo_c2eavgcoeffs( regions(ireg) ) ! cell2edge averaging
308  CALL rflo_checkmetrics( ireg,regions(ireg) ) ! check metrics
309 ! IF (regions(iReg)%levels(1)%grid%remesh==1) THEN
310 ! CALL RFLO_GridRemesh( regions(iReg) ) ! grid remeshing
311 ! iRemesh=1
312 ! ENDIF
313  CALL rflo_calcgridspeeds( regions(ireg) ) ! grid speeds
314  ENDIF ! region on this processor and active, grid moving
315  ENDDO ! iReg
316 
317 #ifdef MPI
318  CALL mpi_allreduce( iremesh, nremesh, 1, mpi_integer, mpi_sum, &
319  global%mpiComm, global%mpierr )
320  IF (global%mpierr /=0 ) CALL errorstop( global,err_mpi_trouble,__line__ )
321  IF (nremesh > 0) someremesh = .true.
322 #endif
323 
324  IF (someremesh) THEN
325  CALL rflo_exchangegeometry( regions ) ! exchange geometry
326  DO ireg=1,global%nRegions
327  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
328  regions(ireg)%active==active .AND. & ! on my processor
329  iremesh==1) THEN ! and remeshing
330  CALL rflo_calcfacevectors( regions(ireg) ) ! faces
331  CALL rflo_calccontrolvolumes( regions(ireg) ) ! volumes
332  CALL rflo_calccellcentroids( regions(ireg) ) ! cell centroids
333  IF (global%moveGridScheme==movegrid_foms) &
334  CALL rflo_calcfacecentroids( regions(ireg) ) ! face centroids
335  IF (regions(ireg)%mixtInput%faceEdgeAvg==fe_avg_linear) &
336  CALL rflo_c2favgcoeffs( regions(ireg) ) ! cell2face averaging
337  CALL rflo_c2eavgcoeffs( regions(ireg) ) ! cell2edge averaging
338  ENDIF ! region on this processor and active, grid moving
339  ENDDO ! iReg
340  ENDIF
341 
342 ! finalize --------------------------------------------------------------------
343 
344  CALL deregisterfunction( global )
345 
346 1000 FORMAT(a,1x,'Global-TFI grid motion:',i6,4(1pe9.2))
347 2000 FORMAT(a,1x,'Global-Weighted-Laplacian grid motion:',i6,4(1pe9.2),1pe13.4)
348 3000 FORMAT(a,1x,'Global-Orthogonal-Laplacian gridmotion:',i4,6(1pe9.2),1pe10.2)
349 4000 FORMAT(a,1x,'global skewness, minvol:',2(1pe14.5))
350 
351 END SUBROUTINE rflo_movegridframe
352 
353 !******************************************************************************
354 !
355 ! Purpose: search for corner points including those of internal patches
356 !
357 ! Description: none.
358 !
359 ! Input: regions = data of current region.
360 !
361 ! Output: grid%nCorns = number of corner points in each region
362 ! grid%ijkCorn = ijkValue of each corner
363 !
364 ! Notes: none.
365 !
366 !******************************************************************************
367 
368 SUBROUTINE rflo_mgframecornpoints( regions )
369 
370  USE modinterfaces, ONLY : rflo_getnodeoffset, &
372 
373  IMPLICIT NONE
374 #include "Indexing.h"
375 
376 ! ... parameters
377  TYPE(t_region), POINTER :: regions(:)
378 
379 ! ... loop variables
380  INTEGER :: l, ipatch, ireg, ipcorn, intcorn, nreg
381 
382 ! ... local variables
383  INTEGER, PARAMETER :: ncmax=100
384  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend
385  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend
386  INTEGER :: iptc, jptc, kptc, iblk, jblk, kblk, ijkcurr
387  INTEGER :: inoff, ijnoff, lbound, regnc, errfl
388  INTEGER, ALLOCATABLE :: ivar(:), ijkcorn(:,:)
389  LOGICAL :: wasfound
390 
391  TYPE(t_patch), POINTER :: patch
392  TYPE(t_grid), POINTER :: grid
393  TYPE(t_global), POINTER :: global
394 
395 !******************************************************************************
396 
397  global => regions(1)%global
398 
399  CALL registerfunction( global,'RFLO_MgFrameCornPoints',&
400  'RFLO_ModMoveGridNconform.F90' )
401 
402 ! search for block and patch corners in each region ---------------------------
403 
404  ALLOCATE( ivar(global%nRegions),stat=errfl )
405  global%error = errfl
406  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
407 
408  ALLOCATE( ijkcorn(ncmax,global%nRegions),stat=errfl )
409  global%error = errfl
410  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
411 
412  ilev = 1
413 
414  DO ireg = 1,global%nRegions
415  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
416  regions(ireg)%active==active) THEN ! on my processor
417 
418  grid => regions(ireg)%levels(ilev)%grid
419 
420  CALL rflo_getdimensphysnodes( regions(ireg),ilev,ipnbeg,ipnend, &
421  jpnbeg,jpnend,kpnbeg,kpnend )
422  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
423 
424 ! --- search for internal patch corners
425 
426  grid%nCorns(ireg) = 8
427  ijkcorn(1,ireg) = indijk(ipnbeg,jpnbeg,kpnbeg,inoff,ijnoff)
428  ijkcorn(2,ireg) = indijk(ipnbeg,jpnbeg,kpnend,inoff,ijnoff)
429  ijkcorn(3,ireg) = indijk(ipnbeg,jpnend,kpnend,inoff,ijnoff)
430  ijkcorn(4,ireg) = indijk(ipnbeg,jpnend,kpnbeg,inoff,ijnoff)
431  ijkcorn(5,ireg) = indijk(ipnend,jpnbeg,kpnbeg,inoff,ijnoff)
432  ijkcorn(6,ireg) = indijk(ipnend,jpnbeg,kpnend,inoff,ijnoff)
433  ijkcorn(7,ireg) = indijk(ipnend,jpnend,kpnend,inoff,ijnoff)
434  ijkcorn(8,ireg) = indijk(ipnend,jpnend,kpnbeg,inoff,ijnoff)
435 
436  DO ipatch=1,regions(ireg)%nPatches
437  patch => regions(ireg)%levels(ilev)%patches(ipatch)
438  lbound = patch%lbound
439 
440  CALL rflo_getpatchindicesnodes( regions(ireg),patch,ilev, &
441  ibeg,iend,jbeg,jend,kbeg,kend )
442 
443  DO ipcorn = 1,4 ! patch corners
444  IF (lbound==1 .OR. lbound==2) THEN
445  iptc = ibeg
446  IF (lbound==1) iblk = ipnbeg
447  IF (lbound==2) iblk = ipnend
448  IF (ipcorn==1) THEN
449  jptc = jbeg
450  jblk = jpnbeg
451  kptc = kbeg
452  kblk = kpnbeg
453  ELSEIF (ipcorn==2) THEN
454  jptc = jbeg
455  jblk = jpnbeg
456  kptc = kend
457  kblk = kpnend
458  ELSEIF (ipcorn==3) THEN
459  jptc = jend
460  jblk = jpnend
461  kptc = kend
462  kblk = kpnend
463  ELSEIF (ipcorn==4) THEN
464  jptc = jend
465  jblk = jpnend
466  kptc = kbeg
467  kblk = kpnbeg
468  ENDIF
469  ELSEIF (lbound==3 .OR. lbound==4) THEN
470  jptc = jbeg
471  IF (lbound==3) jblk = jpnbeg
472  IF (lbound==4) jblk = jpnend
473  IF (ipcorn==1) THEN
474  kptc = kbeg
475  kblk = kpnbeg
476  iptc = ibeg
477  iblk = ipnbeg
478  ELSEIF (ipcorn==2) THEN
479  kptc = kbeg
480  kblk = kpnbeg
481  iptc = iend
482  iblk = ipnend
483  ELSEIF (ipcorn==3) THEN
484  kptc = kend
485  kblk = kpnend
486  iptc = iend
487  iblk = ipnend
488  ELSEIF (ipcorn==4) THEN
489  kptc = kend
490  kblk = kpnend
491  iptc = ibeg
492  iblk = ipnbeg
493  ENDIF ! ipCorn
494  ELSEIF (lbound==5 .OR. lbound==6) THEN
495  kptc = kbeg
496  IF (lbound==5) kblk = kpnbeg
497  IF (lbound==6) kblk = kpnend
498  IF (ipcorn==1) THEN
499  iptc = ibeg
500  iblk = ipnbeg
501  jptc = jbeg
502  jblk = jpnbeg
503  ELSEIF (ipcorn==2) THEN
504  iptc = ibeg
505  iblk = ipnbeg
506  jptc = jend
507  jblk = jpnend
508  ELSEIF (ipcorn==3) THEN
509  iptc = iend
510  iblk = ipnend
511  jptc = jend
512  jblk = jpnend
513  ELSEIF (ipcorn==4) THEN
514  iptc = iend
515  iblk = ipnend
516  jptc = jbeg
517  jblk = jpnbeg
518  ENDIF ! ipCorn
519  ENDIF ! lbound
520 
521  patch%corns(ipcorn) = indijk(iptc,jptc,kptc,inoff,ijnoff)
522 
523  IF (iptc/=iblk .OR. jptc/=jblk .OR. kptc/=kblk) THEN
524  wasfound = .false.
525  ijkcurr = indijk(iptc,jptc,kptc,inoff,ijnoff)
526  DO intcorn=1,grid%nCorns(ireg)
527  IF (ijkcorn(intcorn,ireg)==ijkcurr) THEN
528  wasfound = .true.
529  ENDIF
530  ENDDO
531  IF (.NOT. wasfound) THEN
532  grid%nCorns(ireg) = grid%nCorns(ireg) +1
533  ijkcorn(grid%nCorns(ireg),ireg) = ijkcurr
534  ENDIF
535  ENDIF
536  IF (grid%nCorns(ireg) >= ncmax) THEN
537  CALL errorstop( global,err_illegal_value,__line__, &
538  'too low ncMax in RFLO_ModMoveGridNconform/RFLO_MgFrameCornPoints')
539  ENDIF
540  ENDDO ! ipCorn
541  ENDDO ! iPatch
542 
543  ivar(ireg) = grid%nCorns(ireg)
544  ENDIF ! myProcid
545  ENDDO ! iReg
546 
547 #ifdef MPI
548  DO ireg = 1,global%nRegions
549  CALL mpi_bcast( ivar(ireg),1,mpi_integer, &
550  regions(ireg)%procId,global%mpiComm,global%mpierr )
551  IF (global%mpierr /=0 ) CALL errorstop( global,err_mpi_trouble,__line__ )
552 
553  CALL mpi_bcast( ijkcorn(1:ncmax,ireg),ncmax,mpi_integer, &
554  regions(ireg)%procId,global%mpiComm,global%mpierr )
555  IF (global%mpierr /=0 ) CALL errorstop( global,err_mpi_trouble,__line__ )
556  ENDDO
557 #endif
558 
559  regnc = 0
560  DO ireg = 1,global%nRegions
561  regnc = max( regnc,ivar(ireg) )
562  ENDDO
563  global%moveGridRegNc = regnc
564 
565  DO ireg = 1,global%nRegions
566  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
567  regions(ireg)%active==active) THEN ! on my processor
568 
569  grid => regions(ireg)%levels(ilev)%grid
570  DO nreg = 1,global%nRegions
571  grid%nCorns(nreg) = ivar(nreg)
572  ENDDO
573 
574  ALLOCATE( grid%ijkCorn( regnc,global%nRegions),stat=errfl )
575  global%error = errfl
576  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
577 
578  ALLOCATE( grid%regCorn( 3,regnc,global%nRegions),stat=errfl )
579  global%error = errfl
580  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
581 
582  ALLOCATE( grid%regCornOld( 3,regnc,global%nRegions),stat=errfl )
583  global%error = errfl
584  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
585 
586  ALLOCATE( grid%regCornOrig(3,regnc,global%nRegions),stat=errfl )
587  global%error = errfl
588  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
589 
590  ALLOCATE( grid%nghbor( 3,6,regnc) ,stat=errfl )
591  global%error = errfl
592  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
593 
594  DO l = 1,grid%nCorns(ireg)
595  grid%ijkCorn(l,ireg) = ijkcorn(l,ireg)
596  ENDDO
597 
598  ENDIF ! myProcid
599  ENDDO ! iReg
600 
601 ! deallocate temporary arrays
602 
603  DEALLOCATE( ivar,stat=errfl )
604  global%error = errfl
605  IF (global%error /= 0) CALL errorstop( global,err_deallocate,__line__ )
606 
607 ! finalize --------------------------------------------------------------------
608 
609  CALL deregisterfunction( global )
610 
611 END SUBROUTINE rflo_mgframecornpoints
612 
613 !******************************************************************************
614 !
615 ! Purpose: broadcast movements at 8 corner points of current region to all
616 ! regions
617 !
618 ! Description: none.
619 !
620 ! Input: regions = data of all grid regions.
621 !
622 ! Notes: upon first call by RFLO_InitGridProcedure, regions%levels%grid%xyz
623 ! contains grid coordinates, but on second call by RFLO_MgFrameSurface
624 ! regions%levels%grid%xyz contains grid movements.
625 !
626 !******************************************************************************
627 
628 SUBROUTINE rflo_mgframebroadcast( regions,iselect,iter )
629 
631  IMPLICIT NONE
632 
633 #include "Indexing.h"
634 
635 ! ... parameters
636  TYPE(t_region), POINTER :: regions(:)
637  INTEGER :: iselect, iter
638 
639 ! ... loop variables
640  INTEGER :: i, l, ireg
641 
642 ! ... local variables
643  INTEGER :: ilev, ncorns, errfl
644  INTEGER :: ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend, inoff, ijnoff
645  INTEGER, ALLOCATABLE :: corner(:)
646 
647  REAL(RFREAL), ALLOCATABLE :: rvar(:,:,:)
648  REAL(RFREAL), POINTER :: dxyz(:,:)
649 
650  TYPE(t_global), POINTER :: global
651  TYPE(t_grid), POINTER :: grid, gridold
652 
653 !******************************************************************************
654 
655  global => regions(1)%global
656 
657  CALL registerfunction( global,'RFLO_MgFrameBroadcast',&
658  'RFLO_ModMoveGridNconform.F90' )
659 
660 ! store block corners and broadcast to all regions ----------------------------
661 
662  ALLOCATE( rvar(xcoord:zcoord,global%moveGridRegNc,global%nRegions), &
663  stat=errfl )
664  global%error = errfl
665  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
666  rvar = 0._rfreal
667 
668  ilev = 1
669 
670  DO ireg = 1,global%nRegions
671  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
672  regions(ireg)%active==active) THEN ! on my processor
673 
674  grid => regions(ireg)%levels(ilev)%grid
675 ! gridOld => regions(iReg)%levels(iLev)%gridOld
676 
677  IF (iter==1) THEN
678  CALL rflo_getdimensphysnodes( regions(ireg),ilev,ipnbeg,ipnend, &
679  jpnbeg,jpnend,kpnbeg,kpnend )
680  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
681 
682  dxyz => grid%xyz
683  ncorns = grid%nCorns(ireg)
684 
685  ALLOCATE( corner(ncorns), stat=errfl )
686  global%error = errfl
687  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
688 
689  DO l = 1,ncorns
690  corner(l) = grid%ijkCorn(l,ireg)
691  ENDDO
692 
693  IF (iselect==1) THEN
694  DO i=1,ncorns
695  grid%regCornOld(xcoord,i,ireg) = dxyz(xcoord,corner(i))
696  grid%regCornOld(ycoord,i,ireg) = dxyz(ycoord,corner(i))
697  grid%regCornOld(zcoord,i,ireg) = dxyz(zcoord,corner(i))
698  rvar(:,i,ireg) = grid%regCornOld(:,i,ireg)
699  ENDDO
700  ELSE
701  DO i=1,ncorns
702  grid%regCornOrig(xcoord,i,ireg) = dxyz(xcoord,corner(i))
703  grid%regCornOrig(ycoord,i,ireg) = dxyz(ycoord,corner(i))
704  grid%regCornOrig(zcoord,i,ireg) = dxyz(zcoord,corner(i))
705  rvar(:,i,ireg) = grid%regCornOrig(:,i,ireg)
706  ENDDO
707  ENDIF
708 
709  DEALLOCATE( corner, stat=errfl )
710  global%error = errfl
711  IF (global%error /= 0) CALL errorstop( global,err_deallocate,__line__ )
712 
713  ELSE
714  rvar(:,:,ireg) = grid%regCornOld(:,:,ireg)
715 
716  ENDIF ! iter
717  ENDIF ! myProcid
718  ENDDO ! iReg
719 
720 #ifdef MPI
721  DO ireg = 1,global%nRegions
722  ncorns = global%moveGridRegNc
723 
724  CALL mpi_bcast( rvar(xcoord:zcoord,1:ncorns,ireg),3*ncorns, &
725  mpi_rfreal,regions(ireg)%procId,global%mpiComm,global%mpierr )
726  IF (global%mpierr /=0 ) CALL errorstop( global,err_mpi_trouble,__line__ )
727  ENDDO
728  CALL mpi_barrier( global%mpiComm,global%mpierr )
729 
730  DO ireg = 1,global%nRegions
731  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
732  regions(ireg)%active==active) THEN ! on my processor
733 
734  grid => regions(ireg)%levels(ilev)%grid
735  IF (iter == 1) THEN
736  IF (iselect==1) THEN
737  DO l=1,global%nRegions
738  grid%regCornOld(:,:,l) = rvar(:,:,l)
739  ENDDO
740  ELSE
741  DO l=1,global%nRegions
742  grid%regCornOrig(:,:,l) = rvar(:,:,l)
743  ENDDO
744  ENDIF
745  ELSE
746  DO l=1,global%nRegions
747  grid%regCornOld(:,:,l) = rvar(:,:,l)
748  ENDDO
749  ENDIF ! iter
750  ENDIF ! myProcid
751  ENDDO ! iReg
752 
753 ! DO iReg = 1,global%nRegions
754 ! IF (regions(iReg)%procid==global%myProcid .AND. & ! region active and
755 ! regions(iReg)%active==ACTIVE) THEN ! on my processor
756 ! grid => regions(1)%levels(iLev)%grid
757 ! DO l = 1,global%nRegions
758 ! DO i=1,grid%nCorns(iReg)
759 ! write(*,*)iReg,l,i,grid%regCornOrig(:,i,l)
760 ! ENDDO
761 ! ENDDO
762 ! ENDIF
763 ! ENDDO
764 #endif
765 
766 ! finalize --------------------------------------------------------------------
767 
768  CALL deregisterfunction( global )
769 
770 END SUBROUTINE rflo_mgframebroadcast
771 
772 !******************************************************************************
773 !
774 ! Purpose: search for six closest neighbors
775 !
776 ! Description: none.
777 !
778 ! Input: regions = data of current region.
779 !
780 ! Output: grid%nghbor = neighbouring points identified
781 !
782 ! Notes: none
783 !
784 !******************************************************************************
785 
786 SUBROUTINE rflo_mgframesrchneighbors( regions )
787 
790 
791  IMPLICIT NONE
792 #include "Indexing.h"
793 
794 ! ... parameters
795  TYPE(t_region), POINTER :: regions(:)
796 
797 ! ... loop variables
798  INTEGER :: i, j, k, ipatch, ic, ireg, nc, nreg
799 
800 ! ... local variables
801  INTEGER :: ilev, bctype, ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend
802  INTEGER :: ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend
803  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, ncorns
804  INTEGER :: ijknode(4), ncmin(6), nregmin(6), inoff, ijnoff, lbound, errfl
805  REAL(RFREAL) :: edgelen, ds, tol, distmin(6)
806  REAL(RFREAL), POINTER :: xyz(:,:)
807  REAL(RFREAL), ALLOCATABLE :: dist(:,:)
808 
809  TYPE(t_patch), POINTER :: patch
810  TYPE(t_grid), POINTER :: grid
811  TYPE(t_global), POINTER :: global
812 
813 !******************************************************************************
814 
815  global => regions(1)%global
816 
817  CALL registerfunction( global,'RFLO_MgFrameSrchNeighbors',&
818  'RFLO_ModMoveGridNconform.F90' )
819 
820 ! search for six closest neighbours -------------------------------------------
821 
822  ilev = 1
823 
824  ALLOCATE( dist(global%moveGridRegNc,global%nRegions), stat=errfl )
825  IF (errfl>0) goto 88
826 
827  DO ireg = 1,global%nRegions
828  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
829  regions(ireg)%active==active) THEN ! on my processor
830 
831  grid => regions(ireg)%levels(ilev)%grid
832 
833  ncorns = grid%nCorns(ireg)
834 
835  CALL rflo_getdimensphys( regions(ireg),ilev,ipcbeg,ipcend, &
836  jpcbeg,jpcend,kpcbeg,kpcend )
837  CALL rflo_getdimensphysnodes( regions(ireg),ilev,ipnbeg,ipnend, &
838  jpnbeg,jpnend,kpnbeg,kpnend )
839  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
840 
841  xyz => regions(ireg)%levels(ilev)%grid%xyz
842 
843 ! --- calculate the shortest cell edge
844 
845  edgelen = 1.e+30_rfreal
846 
847  DO k=kpcbeg,kpcend
848  DO j=jpcbeg,jpcend
849  DO i=ipcbeg,ipcend
850  ijknode(1) = indijk(i ,j ,k ,inoff,ijnoff)
851  ijknode(2) = indijk(i+1,j ,k ,inoff,ijnoff)
852  ijknode(3) = indijk(i ,j+1,k ,inoff,ijnoff)
853  ijknode(4) = indijk(i ,j ,k+1,inoff,ijnoff)
854  ds = sqrt((xyz(xcoord,ijknode(2))-xyz(xcoord,ijknode(1)))**2+ &
855  (xyz(ycoord,ijknode(2))-xyz(ycoord,ijknode(1)))**2+ &
856  (xyz(zcoord,ijknode(2))-xyz(zcoord,ijknode(1)))**2)
857  edgelen = min(edgelen,ds)
858  ds = sqrt((xyz(xcoord,ijknode(3))-xyz(xcoord,ijknode(1)))**2+ &
859  (xyz(ycoord,ijknode(3))-xyz(ycoord,ijknode(1)))**2+ &
860  (xyz(zcoord,ijknode(3))-xyz(zcoord,ijknode(1)))**2)
861  edgelen = min(edgelen,ds)
862  ds = sqrt((xyz(xcoord,ijknode(4))-xyz(xcoord,ijknode(1)))**2+ &
863  (xyz(ycoord,ijknode(4))-xyz(ycoord,ijknode(1)))**2+ &
864  (xyz(zcoord,ijknode(4))-xyz(zcoord,ijknode(1)))**2)
865  edgelen = min(edgelen,ds)
866  ENDDO
867  ENDDO
868  ENDDO
869  tol = 1.e-5_rfreal*edgelen
870 
871  DO ic = 1,ncorns
872  distmin(1:6) = 1.e+30_rfreal
873  ncmin(1:6) = 1
874  nregmin(1:6) = 1
875  DO nreg = 1,global%nRegions
876  DO nc = 1,grid%nCorns(nreg)
877  dist(nc,nreg) = sqrt((grid%regCornOrig(xcoord,nc,nreg)- &
878  grid%regCornOrig(xcoord,ic,ireg))**2 + &
879  (grid%regCornOrig(ycoord,nc,nreg)- &
880  grid%regCornOrig(ycoord,ic,ireg))**2 + &
881  (grid%regCornOrig(zcoord,nc,nreg)- &
882  grid%regCornOrig(zcoord,ic,ireg))**2)
883 
884 ! -------- inhibitor check
885 ! IF (dist(nc,nReg)>edgeLen .AND. iReg==12 .AND. ic==5 .AND. &
886 ! (nReg==12 .OR. nReg==14 .OR. nReg==18 .OR. nReg==25 .OR. &
887 ! nReg==26 .OR. nReg==52 .OR. nReg==69)) &
888 ! write(*,*)'i',iReg,ic,nReg,nc,dist(nc,nReg)
889 
890 ! -------- titan4-240blocks check
891 ! IF (dist(nc,nReg)>edgeLen .AND. iReg==120 .AND. ic==6 .AND. &
892 ! (nReg==96 .OR. nReg==37 .OR. nReg==105 .OR. nReg==121 .OR. &
893 ! nReg==117 .OR. nReg==123 .OR. nReg==71 .OR. nReg==72 .OR. &
894 ! nReg==120)) write(*,*)'i',iReg,ic,nReg,nc,dist(nc,nReg)
895 
896  IF (dist(nc,nreg)<distmin(1) .AND. dist(nc,nreg)>edgelen) THEN
897  DO k = 6,2,-1
898  distmin(k) = distmin(k-1)
899  ncmin(k) = ncmin(k-1)
900  nregmin(k) = nregmin(k-1)
901  ENDDO
902  distmin(1) = dist(nc,nreg)
903  ncmin(1) = nc
904  nregmin(1) = nreg
905  ENDIF
906 
907  DO k = 2,6
908  IF (dist(nc,nreg) > (distmin(k-1) + tol) .AND. &
909  dist(nc,nreg) < (distmin(k) - tol)) THEN
910 ! ------------- titan4-240blocks check
911 ! IF (iReg==120 .AND. ic==6) write(*,*)'j',iReg,ic,nReg, &
912 ! nc,k,dist(nc,nReg),distMin(k-1),distMin(k)
913  DO j = 6,k+1,-1
914  distmin(j) = distmin(j-1)
915  ncmin(j) = ncmin(j-1)
916  nregmin(j) = nregmin(j-1)
917  ENDDO
918  distmin(k) = dist(nc,nreg)
919  ncmin(k) = nc
920  nregmin(k) = nreg
921  ENDIF
922  ENDDO
923 
924  ENDDO ! nc
925  ENDDO ! nReg
926 
927  DO k = 1,6
928  nc = ncmin(k)
929  nreg = nregmin(k)
930  grid%nghbor(1,k,ic) = nc
931  grid%nghbor(2,k,ic) = nreg
932  ENDDO ! k
933  ENDDO ! ic
934 
935  grid%nghbor(3,:,:) = 1
936 
937  ENDIF ! myProcid
938  ENDDO ! iReg
939 
940 ! deallocate temporary arrays -------------------------------------------------
941 
942  DEALLOCATE( dist, stat=errfl ); IF (errfl>0) goto 99
943 
944 ! assign internal/external flag to block corners
945 
946  DO ireg = 1,global%nRegions
947  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
948  regions(ireg)%active==active) THEN ! on my processor
949 
950  grid => regions(ireg)%levels(ilev)%grid
951 
952 ! --- corner 1
953  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(1)%interact).OR. &
954  (.NOT. regions(ireg)%levels(ilev)%edgeCells(1)%interact).OR. &
955  (.NOT. regions(ireg)%levels(ilev)%edgeCells(4)%interact).OR. &
956  (.NOT. regions(ireg)%levels(ilev)%edgeCells(9)%interact)) &
957  grid%nghbor(3,1:6,1) = 0
958 
959 ! --- corner 2
960  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(2)%interact).OR. &
961  (.NOT. regions(ireg)%levels(ilev)%edgeCells(1)%interact).OR. &
962  (.NOT. regions(ireg)%levels(ilev)%edgeCells(2)%interact).OR. &
963  (.NOT. regions(ireg)%levels(ilev)%edgeCells(10)%interact)) &
964  grid%nghbor(3,1:6,2) = 0
965 
966 ! --- corner 3
967  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(3)%interact).OR. &
968  (.NOT. regions(ireg)%levels(ilev)%edgeCells(2)%interact).OR. &
969  (.NOT. regions(ireg)%levels(ilev)%edgeCells(3)%interact).OR. &
970  (.NOT. regions(ireg)%levels(ilev)%edgeCells(11)%interact)) &
971  grid%nghbor(3,1:6,3) = 0
972 
973 ! --- corner 4
974  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(4)%interact).OR. &
975  (.NOT. regions(ireg)%levels(ilev)%edgeCells(3)%interact).OR. &
976  (.NOT. regions(ireg)%levels(ilev)%edgeCells(4)%interact).OR. &
977  (.NOT. regions(ireg)%levels(ilev)%edgeCells(12)%interact)) &
978  grid%nghbor(3,1:6,4) = 0
979 
980 ! --- corner 5
981  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(5)%interact).OR. &
982  (.NOT. regions(ireg)%levels(ilev)%edgeCells(5)%interact).OR. &
983  (.NOT. regions(ireg)%levels(ilev)%edgeCells(8)%interact).OR. &
984  (.NOT. regions(ireg)%levels(ilev)%edgeCells(9)%interact)) &
985  grid%nghbor(3,1:6,5) = 0
986 
987 ! --- corner 6
988  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(6)%interact).OR. &
989  (.NOT. regions(ireg)%levels(ilev)%edgeCells(5)%interact).OR. &
990  (.NOT. regions(ireg)%levels(ilev)%edgeCells(6)%interact).OR. &
991  (.NOT. regions(ireg)%levels(ilev)%edgeCells(10)%interact)) &
992  grid%nghbor(3,1:6,6) = 0
993 
994 ! --- corner 7
995  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(7)%interact).OR. &
996  (.NOT. regions(ireg)%levels(ilev)%edgeCells(6)%interact).OR. &
997  (.NOT. regions(ireg)%levels(ilev)%edgeCells(7)%interact).OR. &
998  (.NOT. regions(ireg)%levels(ilev)%edgeCells(11)%interact)) &
999  grid%nghbor(3,1:6,7) = 0
1000 
1001 ! --- corner 8
1002  IF ((.NOT. regions(ireg)%levels(ilev)%cornerCells(8)%interact).OR. &
1003  (.NOT. regions(ireg)%levels(ilev)%edgeCells(7)%interact).OR. &
1004  (.NOT. regions(ireg)%levels(ilev)%edgeCells(8)%interact).OR. &
1005  (.NOT. regions(ireg)%levels(ilev)%edgeCells(12)%interact)) &
1006  grid%nghbor(3,1:6,8) = 0
1007 
1008  ENDIF ! myProcid
1009  ENDDO ! iReg
1010 
1011 ! second step search for external corners
1012 
1013  DO ireg = 1,global%nRegions
1014  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1015  regions(ireg)%active==active) THEN ! on my processor
1016 
1017  grid => regions(ireg)%levels(ilev)%grid
1018 
1019  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
1020 
1021  DO ipatch=1,regions(ireg)%nPatches
1022  patch => regions(ireg)%levels(ilev)%patches(ipatch)
1023  lbound = patch%lbound
1024  bctype = patch%bcType
1025 
1026  CALL rflo_getpatchindicesnodes( regions(ireg),patch,ilev, &
1027  ibeg,iend,jbeg,jend,kbeg,kend )
1028 
1029  IF ((bctype>=bc_inflow .AND. bctype<=bc_inflow +bc_range) .OR. &
1030  (bctype>=bc_outflow .AND. bctype<=bc_outflow +bc_range) .OR. &
1031  (bctype>=bc_slipwall .AND. bctype<=bc_slipwall +bc_range) .OR. &
1032  (bctype>=bc_noslipwall .AND. bctype<=bc_noslipwall+bc_range) .OR. &
1033  (bctype>=bc_farfield .AND. bctype<=bc_farfield +bc_range) .OR. &
1034  (bctype>=bc_injection .AND. bctype<=bc_injection +bc_range) .OR. &
1035  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
1036  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
1037  IF (lbound==1 .OR. lbound==2) THEN
1038  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1039  ijknode(2) = indijk(ibeg,jbeg,kend,inoff,ijnoff)
1040  ijknode(3) = indijk(ibeg,jend,kend,inoff,ijnoff)
1041  ijknode(4) = indijk(ibeg,jend,kbeg,inoff,ijnoff)
1042  ELSEIF (lbound==3 .OR. lbound==4) THEN
1043  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1044  ijknode(2) = indijk(ibeg,jbeg,kend,inoff,ijnoff)
1045  ijknode(3) = indijk(iend,jbeg,kend,inoff,ijnoff)
1046  ijknode(4) = indijk(iend,jbeg,kbeg,inoff,ijnoff)
1047  ELSEIF (lbound==5 .OR. lbound==6) THEN
1048  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1049  ijknode(2) = indijk(ibeg,jend,kbeg,inoff,ijnoff)
1050  ijknode(3) = indijk(iend,jend,kbeg,inoff,ijnoff)
1051  ijknode(4) = indijk(iend,jbeg,kbeg,inoff,ijnoff)
1052  ENDIF ! lbound
1053  DO ic = 1,grid%nCorns(ireg)
1054  IF (ijknode(1)==grid%ijkCorn(ic,ireg) .OR. &
1055  ijknode(2)==grid%ijkCorn(ic,ireg) .OR. &
1056  ijknode(3)==grid%ijkCorn(ic,ireg) .OR. &
1057  ijknode(4)==grid%ijkCorn(ic,ireg)) grid%nghbor(3,1:6,ic)= 2
1058  ENDDO
1059  ENDIF ! bc_external
1060  ENDDO ! iPatch
1061 
1062  ENDIF ! myProcid
1063  ENDDO ! iReg
1064 
1065 ! third step search for external corners
1066 
1067  DO ireg = 1,global%nRegions
1068  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1069  regions(ireg)%active==active) THEN ! on my processor
1070 
1071  grid => regions(ireg)%levels(ilev)%grid
1072 
1073  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
1074 
1075  DO ipatch=1,regions(ireg)%nPatches
1076  patch => regions(ireg)%levels(ilev)%patches(ipatch)
1077  lbound = patch%lbound
1078  bctype = patch%bcType
1079 
1080  CALL rflo_getpatchindicesnodes( regions(ireg),patch,ilev, &
1081  ibeg,iend,jbeg,jend,kbeg,kend )
1082 
1083  IF (patch%bcMotion == bc_external) THEN
1084  IF (lbound==1 .OR. lbound==2) THEN
1085  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1086  ijknode(2) = indijk(ibeg,jbeg,kend,inoff,ijnoff)
1087  ijknode(3) = indijk(ibeg,jend,kend,inoff,ijnoff)
1088  ijknode(4) = indijk(ibeg,jend,kbeg,inoff,ijnoff)
1089  ELSEIF (lbound==3 .OR. lbound==4) THEN
1090  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1091  ijknode(2) = indijk(ibeg,jbeg,kend,inoff,ijnoff)
1092  ijknode(3) = indijk(iend,jbeg,kend,inoff,ijnoff)
1093  ijknode(4) = indijk(iend,jbeg,kbeg,inoff,ijnoff)
1094  ELSEIF (lbound==5 .OR. lbound==6) THEN
1095  ijknode(1) = indijk(ibeg,jbeg,kbeg,inoff,ijnoff)
1096  ijknode(2) = indijk(ibeg,jend,kbeg,inoff,ijnoff)
1097  ijknode(3) = indijk(iend,jend,kbeg,inoff,ijnoff)
1098  ijknode(4) = indijk(iend,jbeg,kbeg,inoff,ijnoff)
1099  ENDIF ! lbound
1100  DO ic = 1,grid%nCorns(ireg)
1101  IF (ijknode(1)==grid%ijkCorn(ic,ireg) .OR. &
1102  ijknode(2)==grid%ijkCorn(ic,ireg) .OR. &
1103  ijknode(3)==grid%ijkCorn(ic,ireg) .OR. &
1104  ijknode(4)==grid%ijkCorn(ic,ireg)) grid%nghbor(3,1:6,ic)= 0
1105  ENDDO
1106  ENDIF ! bc_external
1107  ENDDO ! iPatch
1108 
1109 ! DO ic = 1,grid%nCorns(iReg)
1110 ! DO k=1,6
1111 ! ------- inhibitor
1112 ! IF (iReg==70) &
1113 ! ------- titan4
1114 ! IF (iReg==120) &
1115 ! write(*,*)iReg,ic,k,edgelen,grid%nghbor(1:3,k,ic)
1116 ! ENDDO
1117 ! ENDDO
1118 
1119  ENDIF ! myProcid
1120  ENDDO ! iReg
1121 
1122  goto 999
1123 
1124 ! finalize --------------------------------------------------------------------
1125 
1126 88 CONTINUE
1127 
1128  global%error = errfl
1129  CALL errorstop( global,err_allocate,__line__ )
1130 
1131 99 CONTINUE
1132 
1133  global%error = errfl
1134  CALL errorstop( global,err_deallocate,__line__ )
1135 
1136 999 CONTINUE
1137 
1138  CALL deregisterfunction( global )
1139 
1140 END SUBROUTINE rflo_mgframesrchneighbors
1141 
1142 
1143 !******************************************************************************
1144 !
1145 ! Purpose: correct six closest neighbors
1146 !
1147 ! Description: none.
1148 !
1149 ! Input: regions = data of current region.
1150 !
1151 ! Output: grid%nghbor = corrected neighbouring points identified
1152 !
1153 ! Notes: none
1154 !
1155 !******************************************************************************
1156 
1157 SUBROUTINE rflo_mgframecorrectneighbors( regions )
1158 
1160 
1161  IMPLICIT NONE
1162 #include "Indexing.h"
1163 
1164 ! ... parameters
1165  TYPE(t_region), POINTER :: regions(:)
1166 
1167 ! ... loop variables
1168  INTEGER :: i, j, k, ic, ireg, nc, nreg, lc, lreg
1169 
1170 ! ... local variables
1171  INTEGER :: ilev, bctype, ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend
1172  INTEGER :: ijknode(4), inoff, ijnoff, errfl
1173  REAL(RFREAL) :: edgelen, ds, du2, dumax
1174  REAL(RFREAL), POINTER :: xyz(:,:)
1175  REAL(RFREAL), ALLOCATABLE :: dist(:,:)
1176 
1177  TYPE(t_patch), POINTER :: patch
1178  TYPE(t_grid), POINTER :: grid
1179  TYPE(t_global), POINTER :: global
1180 
1181 !******************************************************************************
1182 
1183  global => regions(1)%global
1184 
1185  CALL registerfunction( global,'RFLO_MgFrameCorrectNeighbors',&
1186  'RFLO_ModMoveGridNconform.F90' )
1187 
1188 ! search for six closest neighbours -------------------------------------------
1189 
1190  ALLOCATE( dist(global%moveGridRegNc,global%nRegions), stat=errfl )
1191  IF (errfl>0) goto 88
1192 
1193  ilev = 1
1194 
1195  DO ireg = 1,global%nRegions
1196  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1197  regions(ireg)%active==active) THEN ! on my processor
1198 
1199  grid => regions(ireg)%levels(ilev)%grid
1200 
1201  CALL rflo_getdimensphys( regions(ireg),ilev,ipcbeg,ipcend, &
1202  jpcbeg,jpcend,kpcbeg,kpcend )
1203  CALL rflo_getnodeoffset( regions(ireg),ilev,inoff,ijnoff )
1204 
1205  xyz => regions(ireg)%levels(ilev)%gridOld%xyz
1206 
1207 ! --- calculate the shortest cell edge
1208 
1209  edgelen = 1.e+30_rfreal
1210 
1211  DO k=kpcbeg,kpcend
1212  DO j=jpcbeg,jpcend
1213  DO i=ipcbeg,ipcend
1214  ijknode(1) = indijk(i ,j ,k ,inoff,ijnoff)
1215  ijknode(2) = indijk(i+1,j ,k ,inoff,ijnoff)
1216  ijknode(3) = indijk(i ,j+1,k ,inoff,ijnoff)
1217  ijknode(4) = indijk(i ,j ,k+1,inoff,ijnoff)
1218  ds = sqrt((xyz(xcoord,ijknode(2))-xyz(xcoord,ijknode(1)))**2+ &
1219  (xyz(ycoord,ijknode(2))-xyz(ycoord,ijknode(1)))**2+ &
1220  (xyz(zcoord,ijknode(2))-xyz(zcoord,ijknode(1)))**2)
1221  edgelen = min(edgelen,ds)
1222  ds = sqrt((xyz(xcoord,ijknode(3))-xyz(xcoord,ijknode(1)))**2+ &
1223  (xyz(ycoord,ijknode(3))-xyz(ycoord,ijknode(1)))**2+ &
1224  (xyz(zcoord,ijknode(3))-xyz(zcoord,ijknode(1)))**2)
1225  edgelen = min(edgelen,ds)
1226  ds = sqrt((xyz(xcoord,ijknode(4))-xyz(xcoord,ijknode(1)))**2+ &
1227  (xyz(ycoord,ijknode(4))-xyz(ycoord,ijknode(1)))**2+ &
1228  (xyz(zcoord,ijknode(4))-xyz(zcoord,ijknode(1)))**2)
1229  edgelen = min(edgelen,ds)
1230  ENDDO
1231  ENDDO
1232  ENDDO
1233 
1234  DO ic = 1,grid%nCorns(ireg)
1235  DO k = 1,6
1236  nc = grid%nghbor(1,k,ic)
1237  nreg = grid%nghbor(2,k,ic)
1238  dumax = -1.e+20_rfreal
1239 
1240  DO lreg = 1,global%nRegions
1241  DO lc = 1,grid%nCorns(lreg)
1242  dist(lc,lreg) = sqrt((grid%regCornOrig(xcoord,nc,nreg)- &
1243  grid%regCornOrig(xcoord,lc,lreg))**2 + &
1244  (grid%regCornOrig(ycoord,nc,nreg)- &
1245  grid%regCornOrig(ycoord,lc,lreg))**2 + &
1246  (grid%regCornOrig(zcoord,nc,nreg)- &
1247  grid%regCornOrig(zcoord,lc,lreg))**2)
1248 
1249  IF (dist(lc,lreg) < 0.1_rfreal*edgelen) THEN
1250  du2 = grid%regCornOld(xcoord,lc,lreg)**2 + &
1251  grid%regCornOld(ycoord,lc,lreg)**2 + &
1252  grid%regCornOld(zcoord,lc,lreg)**2
1253 
1254  IF ( du2 > dumax ) THEN
1255  dumax = du2
1256  grid%nghbor(1,k,ic) = lc
1257  grid%nghbor(2,k,ic) = lreg
1258  ENDIF ! duMax
1259  ENDIF ! dist
1260  ENDDO ! lc
1261  ENDDO ! lReg
1262  ENDDO ! k
1263  ENDDO ! ic
1264 
1265  ENDIF ! myProcid
1266  ENDDO ! iReg
1267 
1268 ! deallocate temporary arrays -------------------------------------------------
1269 
1270  DEALLOCATE( dist, stat=errfl ); IF (errfl>0) goto 99
1271 
1272  goto 999
1273 
1274 ! finalize --------------------------------------------------------------------
1275 
1276 88 CONTINUE
1277 
1278  global%error = errfl
1279  CALL errorstop( global,err_allocate,__line__ )
1280 
1281 99 CONTINUE
1282 
1283  global%error = errfl
1284  CALL errorstop( global,err_deallocate,__line__ )
1285 
1286 999 CONTINUE
1287 
1288  CALL deregisterfunction( global )
1289 
1290 END SUBROUTINE rflo_mgframecorrectneighbors
1291 
1292 !******************************************************************************
1293 !
1294 ! Purpose: move block corners by averaging over 6 closest neighbours
1295 !
1296 ! Description: none.
1297 !
1298 ! Input: regions = data of current region.
1299 !
1300 ! Output: region%levels%grid%regCorn = new block corners movement.
1301 !
1302 ! Notes: none
1303 !
1304 !******************************************************************************
1305 
1306 SUBROUTINE rflo_mgframemovecorners( regions )
1307 
1308  USE modtools, ONLY : isnan
1309  IMPLICIT NONE
1310 
1311 #include "Indexing.h"
1312 
1313 ! ... parameters
1314  TYPE(t_region), POINTER :: regions(:)
1315 
1316 ! ... loop variables
1317  INTEGER :: ireg, ico, k
1318 
1319 ! ... local variables
1320  INTEGER :: ilev, interior, nco(6), nreg(6), ijkcorn
1321  INTEGER :: ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend, inoff, ijnoff
1322  REAL(RFREAL) :: rdenom, amp(3), pow, dist(6), wght(6)
1323 
1324  TYPE(t_grid), POINTER :: grid
1325  TYPE(t_global), POINTER :: global
1326 
1327 !******************************************************************************
1328 
1329  global => regions(1)%global
1330 
1331  CALL registerfunction( global,'RFLO_MgFrameMoveCorners',&
1332  'RFLO_ModMoveGridNconform.F90' )
1333 
1334 ! move block corners ----------------------------------------------------------
1335 
1336  ilev = 1
1337  amp(1) = global%moveGridAmplifX
1338  amp(2) = global%moveGridAmplifY
1339  amp(3) = global%moveGridAmplifZ
1340  pow = global%moveGridPower
1341 
1342  DO ireg = 1,global%nRegions
1343  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1344  regions(ireg)%active==active) THEN ! on my processor
1345 
1346  grid => regions(ireg)%levels(ilev)%grid
1347 
1348  DO ico = 1,grid%nCorns(ireg)
1349  nco(1:6) = grid%nghbor(1,1:6,ico)
1350  nreg(1:6) = grid%nghbor(2,1:6,ico)
1351  interior = grid%nghbor(3,1 ,ico)
1352 
1353  IF (interior==1) THEN
1354  DO k = 1,6
1355  dist(k) = (grid%regCornOrig(xcoord,nco(k),nreg(k)) - &
1356  grid%regCornOrig(xcoord,ico,ireg))**2 + &
1357  (grid%regCornOrig(ycoord,nco(k),nreg(k)) - &
1358  grid%regCornOrig(ycoord,ico,ireg))**2 + &
1359  (grid%regCornOrig(zcoord,nco(k),nreg(k)) - &
1360  grid%regCornOrig(zcoord,ico,ireg))**2
1361  dist(k) = 1._rfreal/sqrt( dist(k) )**pow
1362  ENDDO
1363  rdenom = 1._rfreal/(dist(1)+dist(2)+dist(3)+dist(4)+dist(5)+dist(6))
1364 
1365  DO k = 1,6
1366  wght(k) = dist(k)*rdenom
1367 ! write(*,*)iReg,ico,k,nReg(k),nco(k),dist(k),wght(k)
1368  IF (isnan(wght(k))) &
1369  CALL errorstop( global,err_illegal_value,__line__, &
1370  'invalid weights for global frame motion')
1371  ENDDO
1372 
1373  grid%regCorn(xcoord,ico,ireg) = &
1374  (wght(1)*grid%regCornOld(xcoord,nco(1),nreg(1)) + &
1375  wght(2)*grid%regCornOld(xcoord,nco(2),nreg(2)) + &
1376  wght(3)*grid%regCornOld(xcoord,nco(3),nreg(3)) + &
1377  wght(4)*grid%regCornOld(xcoord,nco(4),nreg(4)) + &
1378  wght(5)*grid%regCornOld(xcoord,nco(5),nreg(5)) + &
1379  wght(6)*grid%regCornOld(xcoord,nco(6),nreg(6)))
1380 
1381  grid%regCorn(ycoord,ico,ireg) = &
1382  (wght(1)*grid%regCornOld(ycoord,nco(1),nreg(1)) + &
1383  wght(2)*grid%regCornOld(ycoord,nco(2),nreg(2)) + &
1384  wght(3)*grid%regCornOld(ycoord,nco(3),nreg(3)) + &
1385  wght(4)*grid%regCornOld(ycoord,nco(4),nreg(4)) + &
1386  wght(5)*grid%regCornOld(ycoord,nco(5),nreg(5)) + &
1387  wght(6)*grid%regCornOld(ycoord,nco(6),nreg(6)))
1388 
1389  grid%regCorn(zcoord,ico,ireg) = &
1390  (wght(1)*grid%regCornOld(zcoord,nco(1),nreg(1)) + &
1391  wght(2)*grid%regCornOld(zcoord,nco(2),nreg(2)) + &
1392  wght(3)*grid%regCornOld(zcoord,nco(3),nreg(3)) + &
1393  wght(4)*grid%regCornOld(zcoord,nco(4),nreg(4)) + &
1394  wght(5)*grid%regCornOld(zcoord,nco(5),nreg(5)) + &
1395  wght(6)*grid%regCornOld(zcoord,nco(6),nreg(6)))
1396  ENDIF ! interior
1397  ENDDO ! ico
1398  ENDIF ! myProcid
1399  ENDDO ! iReg
1400 
1401  DO ireg = 1,global%nRegions
1402  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1403  regions(ireg)%active==active) THEN ! on my processor
1404 
1405  grid => regions(ireg)%levels(ilev)%grid
1406 
1407  DO ico = 1,grid%nCorns(ireg)
1408  interior = grid%nghbor(3, 1, ico)
1409  IF (interior==1) THEN
1410  grid%regCornOld(xcoord,ico,ireg)=amp(1)*grid%regCorn(xcoord,ico,ireg)
1411  grid%regCornOld(ycoord,ico,ireg)=amp(2)*grid%regCorn(ycoord,ico,ireg)
1412  grid%regCornOld(zcoord,ico,ireg)=amp(3)*grid%regCorn(zcoord,ico,ireg)
1413 
1414  ijkcorn = grid%ijkCorn(ico,ireg)
1415  grid%xyz(xcoord,ijkcorn) = grid%regCorn(xcoord,ico,ireg)
1416  grid%xyz(ycoord,ijkcorn) = grid%regCorn(ycoord,ico,ireg)
1417  grid%xyz(zcoord,ijkcorn) = grid%regCorn(zcoord,ico,ireg)
1418  ENDIF
1419  ENDDO ! ico
1420 
1421  ENDIF ! myProcid
1422  ENDDO ! iReg
1423 
1424 ! finalize --------------------------------------------------------------------
1425 
1426  CALL deregisterfunction( global )
1427 
1428 END SUBROUTINE rflo_mgframemovecorners
1429 
1430 
1431 !******************************************************************************
1432 !
1433 ! Purpose: receive and distribute the deformations of surfaces
1434 ! in block-wise manner.
1435 !
1436 ! Description: none.
1437 !
1438 ! Input: regions = data of all grid regions.
1439 !
1440 ! Output: regions%levels%grid%xyz = deformations at the boundaries
1441 ! someMoved = parts of grid moved.
1442 !
1443 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
1444 ! is applied to the finest grid first.
1445 !
1446 !******************************************************************************
1447 
1448 SUBROUTINE rflo_mgframesurfaces( regions,someMoved,iType )
1449 
1453  IMPLICIT NONE
1454 
1455 ! ... parameters
1456  LOGICAL :: somemoved
1457  INTEGER :: itype
1458 
1459  TYPE(t_region), POINTER :: regions(:)
1460 
1461 ! ... loop variables
1462  INTEGER :: ireg, iter, ipatch, i, j, k, ijkn
1463 
1464 ! ... local variables
1465  INTEGER :: ilev, bctype
1466  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, inoff, ijnoff
1467  TYPE(t_grid), POINTER :: grid, gridold
1468  TYPE(t_global), POINTER :: global
1469  TYPE(t_patch), POINTER :: patch
1470 
1471 !******************************************************************************
1472 
1473  global => regions(1)%global
1474 
1475  CALL registerfunction( global,'RFLO_MgFrameSurfaces',&
1476  'RFLO_ModMoveGridNconform.F90' )
1477 
1478 ! move grid separately for each region ----------------------------------------
1479 
1480  somemoved = .false.
1481  ilev = 1
1482 
1483  DO ireg=1,global%nRegions
1484  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1485  regions(ireg)%active==active .AND. & ! on my processor
1486  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
1487 
1488  grid => regions(ireg)%levels(ilev)%grid
1489  gridold => regions(ireg)%levels(ilev)%gridOld
1490  somemoved = .true.
1491 
1492 ! --- store the old grid
1493 
1494  gridold%indSvel = grid%indSvel
1495  gridold%ipc = grid%ipc
1496  gridold%jpc = grid%jpc
1497  gridold%kpc = grid%kpc
1498  gridold%xyz(:,:) = grid%xyz(:,:)
1499  gridold%si(:,:) = grid%si(:,:)
1500  gridold%sj(:,:) = grid%sj(:,:)
1501  gridold%sk(:,:) = grid%sk(:,:)
1502  gridold%vol(:) = grid%vol(:)
1503 
1504 ! --- calculate arclengths between boundaries
1505 
1506  CALL rflo_arclengthbounds( regions(ireg),gridold%xyzOld, &
1507  grid%arcLen12,grid%arcLen34,grid%arcLen56 )
1508 
1509 ! --- get the boundary deformations
1510 
1511  CALL rflo_getdeformation( regions(ireg),grid%boundMoved,grid%xyz )
1512 
1513  grid%xyzOld(:,:) = grid%xyz(:,:)
1514 
1515  ENDIF ! region on this processor and active, grid moving
1516  ENDDO ! iReg
1517 
1518 ! broadcast and compute block corners deformation
1519 
1520  iter = 1
1521  CALL rflo_mgframebroadcast( regions,1,iter )
1522  CALL rflo_mgframecorrectneighbors( regions )
1523  CALL rflo_mgframemovecorners( regions )
1524 
1525  DO iter = 2,10
1526  CALL rflo_mgframebroadcast( regions,1,iter )
1527  CALL rflo_mgframemovecorners( regions )
1528  ENDDO
1529 
1530  DO ireg=1,global%nRegions
1531  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1532  regions(ireg)%active==active .AND. & ! on my processor
1533  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
1534 
1535  grid => regions(ireg)%levels(ilev)%grid
1536  gridold => regions(ireg)%levels(ilev)%gridOld
1537 
1538 ! --- calculate deformations at remaining edges
1539 
1540  CALL rflo_mgframeedgeso( regions(ireg),itype,grid%boundMoved, &
1541  grid%allExternal,grid%edgeMoved,grid%arcLen12, &
1542  grid%arcLen34,grid%arcLen56,gridold%xyzOld,grid%xyz )
1543 ! CALL RFLO_MgFrameEdges( regions(iReg),grid%edgeMoved, &
1544 ! gridOld%xyzOld,grid%xyz )
1545 
1546 ! --- calculate deformations at remaining boundaries
1547 
1548  IF (itype==1) THEN
1549  CALL rflo_mgframebnddeformation( regions(ireg),grid%boundMoved, &
1550  grid%edgeMoved,grid%arcLen12, &
1551  grid%arcLen34,grid%arcLen56, &
1552  gridold%xyzOld,grid%xyz )
1553  ELSE ! iType
1554  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
1555  grid%edgeMoved,grid%arcLen12, &
1556  grid%arcLen34,grid%arcLen56, &
1557  gridold%xyzOld,grid%xyz )
1558  ENDIF ! iType
1559  CALL rflo_mgframerestoreexternal( regions(ireg) )
1560 
1561  ENDIF ! region on this processor and active, grid moving
1562  ENDDO ! iReg
1563 
1564 ! finalize --------------------------------------------------------------------
1565 
1566  CALL deregisterfunction( global )
1567 
1568 END SUBROUTINE rflo_mgframesurfaces
1569 
1570 
1571 !******************************************************************************
1572 !
1573 ! Purpose: restore deformation of solid surfaces from genx at given patches.
1574 !
1575 ! Description: none.
1576 !
1577 ! Input: region = data of current region.
1578 !
1579 ! Output: regions%levels%grid%xyz = deformations at the boundaries restored
1580 !
1581 ! Notes: grid%xyz temporarily stores nodal displacements. The 'untouched'
1582 ! deformation from genx has been saved in grid%xyzOld.
1583 !
1584 !******************************************************************************
1585 
1586 SUBROUTINE rflo_mgframerestoreexternal( region )
1587 
1589  IMPLICIT NONE
1590 #include "Indexing.h"
1591 
1592 ! ... parameters
1593  TYPE(t_region) :: region
1594 
1595 ! ... loop variables
1596  INTEGER :: ireg, ipatch, i, j, k
1597 
1598 ! ... local variables
1599  INTEGER :: ilev, ijkn, lbound
1600  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, inoff, ijnoff
1601  TYPE(t_grid), POINTER :: grid
1602  TYPE(t_global), POINTER :: global
1603  TYPE(t_patch), POINTER :: patch
1604 
1605 !******************************************************************************
1606 
1607  global => region%global
1608 
1609  CALL registerfunction( global,'RFLO_MgFrameRestoreExternal',&
1610  'RFLO_ModMoveGridNconform.F90' )
1611 
1612 ! parameters and pointers -----------------------------------------------------
1613 
1614  ilev = 1
1615  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
1616 
1617  grid => region%levels(ilev)%grid
1618 
1619 ! restore displacements
1620 
1621  DO ipatch=1,region%nPatches
1622  patch => region%levels(ilev)%patches(ipatch)
1623  lbound = patch%lbound
1624 
1625  IF (patch%bcMotion == bc_external .AND. &
1626  (grid%allExternal(lbound).EQV..false.)) THEN
1627 ! IF (patch%bcMotion == BC_EXTERNAL) THEN
1628  CALL rflo_getpatchindicesnodes( region,patch,ilev, &
1629  ibeg,iend,jbeg,jend,kbeg,kend )
1630 
1631  DO k=kbeg,kend
1632  DO j=jbeg,jend
1633  DO i=ibeg,iend
1634  ijkn = indijk(i,j,k,inoff,ijnoff)
1635  grid%xyz(xcoord,ijkn) = grid%xyzOld(xcoord,ijkn)
1636  grid%xyz(ycoord,ijkn) = grid%xyzOld(ycoord,ijkn)
1637  grid%xyz(zcoord,ijkn) = grid%xyzOld(zcoord,ijkn)
1638  ENDDO
1639  ENDDO
1640  ENDDO
1641 
1642  ENDIF ! external BC
1643  ENDDO ! iPatch
1644 
1645 ! finalize --------------------------------------------------------------------
1646 
1647  CALL deregisterfunction( global )
1648 
1649 END SUBROUTINE rflo_mgframerestoreexternal
1650 
1651 !******************************************************************************
1652 !
1653 ! Purpose: calculate node displacements on those edges whose end points have
1654 ! moved, but the associated boundaries were not updated yet (finest
1655 ! grid only).
1656 !
1657 ! Description: points along an edge are shifted using 1-D linear transfinite
1658 ! interpolation (TFI).
1659 !
1660 ! Input: region = grid dimensions
1661 ! boundMoved = flag for boundaries of a region which have moved
1662 ! arcLen12 = arclength between i=const. boundaries for each j, k
1663 ! arcLen34 = arclength between j=const. boundaries for each k, i
1664 ! arcLen56 = arclength between k=const. boundaries for each i, j
1665 ! xyzOld = grid from previous time step.
1666 !
1667 ! Output: edgeMoved = flag if discretization at an edge was changed
1668 ! dNode = updated deformations at edges.
1669 !
1670 ! Notes: variable dNode contains the whole 3-D field.
1671 !
1672 !******************************************************************************
1673 
1674 SUBROUTINE rflo_mgframeedgeso( region,iType,boundMoved,allExternal,edgeMoved, &
1675  arclen12,arclen34,arclen56,xyzold,dnode )
1676 
1678  rflo_tfint1d
1679  IMPLICIT NONE
1680 
1681 #include "Indexing.h"
1682 
1683 ! ... parameters
1684  LOGICAL :: boundmoved(6), allexternal(6), edgemoved(12)
1685 
1686  INTEGER :: itype
1687  REAL(RFREAL), POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
1688  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:)
1689 
1690  TYPE(t_region) :: region
1691 
1692 ! ... loop variables
1693  INTEGER :: iedge, ind
1694 
1695 ! ... local variables
1696  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend, l1c, l2c
1697  INTEGER :: indbeg, indend, ijkn, ijkn1, ijknbeg, ijknend, inoff, ijnoff
1698  INTEGER :: switch(12,11), intertype, iedgeglo
1699 
1700  REAL(RFREAL) :: arclen, ds, s, dn(3), dnbeg(3), dnend(3)
1701  LOGICAL :: interact
1702 
1703 !******************************************************************************
1704 
1705  CALL registerfunction( region%global,'RFLO_MgFrameEdgesO',&
1706  'RFLO_ModMoveGridNconform.F90' )
1707 
1708 ! get dimensions --------------------------------------------------------------
1709 
1710  ilev = 1
1711  CALL rflo_getdimensphysnodes( region,ilev,ipnbeg,ipnend, &
1712  jpnbeg,jpnend,kpnbeg,kpnend )
1713  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
1714 
1715 ! set edge switch -------------------------------------------------------------
1716 ! switch(:,1) = begins at boundary
1717 ! switch(:,2) = ends on boundary
1718 ! switch(:,3) = right boundary
1719 ! switch(:,4) = left boundary
1720 ! switch(:,5) = direction (from-to boundary)
1721 ! switch(:,6) = start index
1722 ! switch(:,7) = end index
1723 ! switch(:,8) = constant index in 1st direction
1724 ! switch(:,9) = constant index in 2nd direction
1725 ! switch(:,10) = start corner number
1726 ! switch(:,11) = end corner number
1727 
1728  switch( 1,:) = (/5, 6, 1, 3, 56, kpnbeg, kpnend, ipnbeg, jpnbeg, 1, 2/)
1729  switch( 2,:) = (/3, 4, 1, 6, 34, jpnbeg, jpnend, kpnend, ipnbeg, 2, 3/)
1730  switch( 3,:) = (/5, 6, 1, 4, 56, kpnbeg, kpnend, ipnbeg, jpnend, 4, 3/)
1731  switch( 4,:) = (/3, 4, 1, 5, 34, jpnbeg, jpnend, kpnbeg, ipnbeg, 1, 4/)
1732  switch( 5,:) = (/5, 6, 2, 3, 56, kpnbeg, kpnend, ipnend, jpnbeg, 5, 6/)
1733  switch( 6,:) = (/3, 4, 2, 6, 34, jpnbeg, jpnend, kpnend, ipnend, 6, 7/)
1734  switch( 7,:) = (/5, 6, 2, 4, 56, kpnbeg, kpnend, ipnend, jpnend, 8, 7/)
1735  switch( 8,:) = (/3, 4, 2, 5, 34, jpnbeg, jpnend, kpnbeg, ipnend, 5, 8/)
1736  switch( 9,:) = (/1, 2, 3, 5, 12, ipnbeg, ipnend, jpnbeg, kpnbeg, 1, 5/)
1737  switch(10,:) = (/1, 2, 3, 6, 12, ipnbeg, ipnend, jpnbeg, kpnend, 2, 6/)
1738  switch(11,:) = (/1, 2, 4, 5, 12, ipnbeg, ipnend, jpnend, kpnbeg, 4, 8/)
1739  switch(12,:) = (/1, 2, 4, 6, 12, ipnbeg, ipnend, jpnend, kpnend, 3, 7/)
1740 
1741 ! edge movement flag ----------------------------------------------------------
1742 
1743  edgemoved(:) = .false.
1744 
1745  IF (itype/=1) THEN
1746  IF (boundmoved(1) .AND. allexternal(1)) THEN
1747  edgemoved( 1) = .true.; edgemoved( 2) = .true.
1748  edgemoved( 3) = .true.; edgemoved( 4) = .true.
1749  ENDIF
1750  IF (boundmoved(2) .AND. allexternal(2)) THEN
1751  edgemoved( 5) = .true.; edgemoved( 6) = .true.
1752  edgemoved( 7) = .true.; edgemoved( 8) = .true.
1753  ENDIF
1754  IF (boundmoved(3) .AND. allexternal(3)) THEN
1755  edgemoved( 1) = .true.; edgemoved( 5) = .true.
1756  edgemoved( 9) = .true.; edgemoved(10) = .true.
1757  ENDIF
1758  IF (boundmoved(4) .AND. allexternal(4)) THEN
1759  edgemoved( 3) = .true.; edgemoved( 7) = .true.
1760  edgemoved(11) = .true.; edgemoved(12) = .true.
1761  ENDIF
1762  IF (boundmoved(5) .AND. allexternal(5)) THEN
1763  edgemoved( 4) = .true.; edgemoved( 8) = .true.
1764  edgemoved( 9) = .true.; edgemoved(11) = .true.
1765  ENDIF
1766  IF (boundmoved(6) .AND. allexternal(6)) THEN
1767  edgemoved( 2) = .true.; edgemoved( 6) = .true.
1768  edgemoved(10) = .true.; edgemoved(12) = .true.
1769  ENDIF
1770  ENDIF ! iType
1771 
1772 ! loop over all 12 edges ------------------------------------------------------
1773 
1774  DO iedge=1,12
1775  IF (.NOT.edgemoved(iedge)) THEN
1776 
1777  edgemoved(iedge) = .true.
1778 
1779  ds = 0._rfreal
1780  indbeg = switch(iedge,6)
1781  indend = switch(iedge,7)
1782  l1c = switch(iedge,8)
1783  l2c = switch(iedge,9)
1784 
1785  iedgeglo = iedge
1786  IF (iedge==11) iedgeglo=12
1787  IF (iedge==12) iedgeglo=11
1788  interact = region%levels(ilev)%edgeCells(iedgeglo)%interact
1789  intertype = region%levels(ilev)%edgeCells(iedgeglo)%interType
1790 
1791  IF (((region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==1 .OR. &
1792  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==1) .AND. &
1793  ((interact .EQV. .true.) .AND. (intertype==edge_interact_full))) &
1794  .OR. &
1795  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==2 .OR. &
1796  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==2) THEN
1797 
1798 ! IF (region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,10))==1 .OR. &
1799 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,11))==1 .OR. &
1800 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,10))==2 .OR. &
1801 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,11))==2) THEN
1802 
1803  DO ind=indbeg+1,indend-1
1804  IF (switch(iedge,5) == 12) THEN
1805  ijkn = indijk(ind ,l1c,l2c,inoff,ijnoff)
1806  ijkn1 = indijk(ind-1 ,l1c,l2c,inoff,ijnoff)
1807  ijknbeg = indijk(indbeg,l1c,l2c,inoff,ijnoff)
1808  ijknend = indijk(indend,l1c,l2c,inoff,ijnoff)
1809  arclen = arclen12(l1c,l2c)
1810  dnbeg(:) = dnode(:,ijknbeg)
1811  dnend(:) = dnode(:,ijknend)
1812  ELSE IF (switch(iedge,5) == 34) THEN
1813  ijkn = indijk(l2c,ind ,l1c,inoff,ijnoff)
1814  ijkn1 = indijk(l2c,ind-1 ,l1c,inoff,ijnoff)
1815  ijknbeg = indijk(l2c,indbeg,l1c,inoff,ijnoff)
1816  ijknend = indijk(l2c,indend,l1c,inoff,ijnoff)
1817  arclen = arclen34(l1c,l2c)
1818  dnbeg(:) = dnode(:,ijknbeg)
1819  dnend(:) = dnode(:,ijknend)
1820  ELSE IF (switch(iedge,5) == 56) THEN
1821  ijkn = indijk(l1c,l2c,ind ,inoff,ijnoff)
1822  ijkn1 = indijk(l1c,l2c,ind-1 ,inoff,ijnoff)
1823  ijknbeg = indijk(l1c,l2c,indbeg,inoff,ijnoff)
1824  ijknend = indijk(l1c,l2c,indend,inoff,ijnoff)
1825  arclen = arclen56(l1c,l2c)
1826  dnbeg(:) = dnode(:,ijknbeg)
1827  dnend(:) = dnode(:,ijknend)
1828  ENDIF
1829  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
1830  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
1831  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
1832  s = ds/arclen
1833 
1834  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
1835  dnode(:,ijkn) = dn(:)
1836  ENDDO ! i
1837  ENDIF ! nghbor
1838  ENDIF ! edgeMoved
1839  ENDDO ! iEdge
1840 
1841 ! finalize --------------------------------------------------------------------
1842 
1843  CALL deregisterfunction( region%global )
1844 
1845 END SUBROUTINE rflo_mgframeedgeso
1846 
1847 !******************************************************************************
1848 !
1849 ! Purpose: calculate node displacements on those patch edges whose end points
1850 ! both are or one of them is interior point.
1851 !
1852 ! Description: points along patch edge are shifted using 1-D linear transfinite
1853 ! interpolation (TFI).
1854 !
1855 ! Input: region = grid dimensions
1856 ! boundMoved = flag for boundaries of a region which have moved
1857 ! xyzOld = grid from previous time step.
1858 !
1859 ! Output: edgeMoved = flag if discretization at an edge was changed
1860 ! dNode = updated deformations at edges.
1861 !
1862 ! Notes: variable dNode contains the whole 3-D field.
1863 !
1864 !******************************************************************************
1865 
1866 SUBROUTINE rflo_mgframeedges( region,edgeMoved,xyzOld,dNode )
1867 
1869  rflo_tfint1d
1870  IMPLICIT NONE
1871 
1872 #include "Indexing.h"
1873 
1874 ! ... parameters
1875  TYPE(t_region) :: region
1876  LOGICAL :: edgemoved(12)
1877  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:)
1878 
1879 ! ... loop variables
1880  INTEGER :: ipedge, ind, ipatch, ic
1881 
1882 ! ... local variables
1883  INTEGER :: ilev, ibeg, iend, jbeg, jend, kbeg, kend, ib, ie, jb, je , kb, ke
1884  INTEGER :: indbeg, indend, ijkn, ijkn1, ijknbeg, ijknend, inoff, ijnoff
1885  INTEGER :: lbound, intb, inte, iedge, intertype
1886 
1887  REAL(RFREAL) :: arclen, ds, s, dn(3), dnbeg(3), dnend(3)
1888  TYPE(t_patch), POINTER :: patch
1889  TYPE(t_grid), POINTER :: grid
1890 
1891 !******************************************************************************
1892 
1893  CALL registerfunction( region%global,'RFLO_MgFrameEdges',&
1894  'RFLO_ModMoveGridNconform.F90' )
1895 
1896 ! get dimensions and pointers -------------------------------------------------
1897 
1898  ilev = 1
1899  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
1900 
1901  grid => region%levels(ilev)%grid
1902 
1903 ! edge movement flag ----------------------------------------------------------
1904 
1905  edgemoved(:) = .true.
1906 
1907 ! loop over all patch edges ---------------------------------------------------
1908 
1909  DO ipatch=1,region%nPatches
1910  patch => region%levels(ilev)%patches(ipatch)
1911  lbound = patch%lbound
1912 
1913  CALL rflo_getpatchindicesnodes( region,patch,ilev, &
1914  ibeg,iend,jbeg,jend,kbeg,kend )
1915 
1916  IF (patch%bcMotion == bc_external) goto 777
1917  DO ipedge = 1,4
1918 
1919  IF (lbound==1 .OR. lbound==2) THEN
1920  IF (ipedge==1) THEN
1921  ib = ibeg
1922  ie = ibeg
1923  jb = jbeg
1924  je = jbeg
1925  kb = kbeg
1926  ke = kend
1927  indbeg = kb
1928  indend = ke
1929  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
1930  ijknend = indijk(ie,je,ke,inoff,ijnoff)
1931  intb = 1
1932  inte = 1
1933  IF (lbound==1) THEN
1934  iedge = 1
1935  ELSE
1936  iedge = 5
1937  ENDIF
1938  intertype = region%levels(ilev)%edgeCells(iedge)%interType
1939  DO ic = 1,grid%nCorns(region%iRegionGlobal)
1940  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
1941  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
1942  ENDDO
1943 
1944  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
1945  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
1946  (intb==2 .OR. inte==2)) THEN
1947 
1948  dnbeg(:) = dnode(:,ijknbeg)
1949  dnend(:) = dnode(:,ijknend)
1950 
1951  arclen = 0._rfreal
1952  DO ind=indbeg+1,indend
1953  ijkn = indijk(ibeg,jbeg,ind ,inoff,ijnoff)
1954  ijkn1 = indijk(ibeg,jbeg,ind-1,inoff,ijnoff)
1955  arclen = arclen + &
1956  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
1957  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
1958  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
1959  ENDDO
1960  ds = 0._rfreal
1961  DO ind=indbeg+1,indend-1
1962  ijkn = indijk(ibeg,jbeg,ind ,inoff,ijnoff)
1963  ijkn1 = indijk(ibeg,jbeg,ind-1,inoff,ijnoff)
1964  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
1965  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
1966  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
1967  s = ds/arclen
1968  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
1969  dnode(:,ijkn) = dn(:)
1970  ENDDO
1971  ENDIF
1972  ELSEIF (ipedge==2) THEN
1973  ib = ibeg
1974  ie = ibeg
1975  jb = jbeg
1976  je = jend
1977  kb = kend
1978  ke = kend
1979  indbeg = jb
1980  indend = je
1981  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
1982  ijknend = indijk(ie,je,ke,inoff,ijnoff)
1983  intb = 1
1984  inte = 1
1985  IF (lbound==1) THEN
1986  iedge = 2
1987  ELSE
1988  iedge = 6
1989  ENDIF
1990  intertype = region%levels(ilev)%edgeCells(iedge)%interType
1991  DO ic = 1,grid%nCorns(region%iRegionGlobal)
1992  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
1993  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
1994  ENDDO
1995 
1996  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
1997  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
1998  (intb==2 .OR. inte==2)) THEN
1999 
2000  dnbeg(:) = dnode(:,ijknbeg)
2001  dnend(:) = dnode(:,ijknend)
2002 
2003  arclen = 0._rfreal
2004  DO ind=indbeg+1,indend
2005  ijkn = indijk(ibeg,ind ,kend,inoff,ijnoff)
2006  ijkn1 = indijk(ibeg,ind-1,kend,inoff,ijnoff)
2007  arclen = arclen + &
2008  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2009  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2010  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2011  ENDDO
2012  ds = 0._rfreal
2013  DO ind=indbeg+1,indend-1
2014  ijkn = indijk(ibeg,ind ,kend,inoff,ijnoff)
2015  ijkn1 = indijk(ibeg,ind-1,kend,inoff,ijnoff)
2016  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2017  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2018  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2019  s = ds/arclen
2020  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2021  dnode(:,ijkn) = dn(:)
2022  ENDDO
2023  ENDIF
2024  ELSEIF (ipedge==3) THEN
2025  ib = ibeg
2026  ie = ibeg
2027  jb = jend
2028  je = jend
2029  kb = kbeg
2030  ke = kend
2031  indbeg = kb
2032  indend = ke
2033  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2034  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2035  intb = 1
2036  inte = 1
2037  IF (lbound==1) THEN
2038  iedge = 3
2039  ELSE
2040  iedge = 7
2041  ENDIF
2042  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2043  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2044  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2045  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2046  ENDDO
2047 
2048  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2049  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2050  (intb==2 .OR. inte==2)) THEN
2051 
2052  dnbeg(:) = dnode(:,ijknbeg)
2053  dnend(:) = dnode(:,ijknend)
2054 
2055  arclen = 0._rfreal
2056  DO ind=indbeg+1,indend
2057  ijkn = indijk(ibeg,jend,ind ,inoff,ijnoff)
2058  ijkn1 = indijk(ibeg,jend,ind-1,inoff,ijnoff)
2059  arclen = arclen + &
2060  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2061  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2062  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2063  ENDDO
2064  ds = 0._rfreal
2065  DO ind=indbeg+1,indend-1
2066  ijkn = indijk(ibeg,jend,ind ,inoff,ijnoff)
2067  ijkn1 = indijk(ibeg,jend,ind-1,inoff,ijnoff)
2068  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2069  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2070  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2071  s = ds/arclen
2072  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2073  dnode(:,ijkn) = dn(:)
2074  ENDDO
2075  ENDIF
2076  ELSEIF (ipedge==4) THEN
2077  ib = ibeg
2078  ie = ibeg
2079  jb = jbeg
2080  je = jend
2081  kb = kbeg
2082  ke = kbeg
2083  indbeg = jb
2084  indend = je
2085  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2086  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2087  intb = 1
2088  inte = 1
2089  IF (lbound==1) THEN
2090  iedge = 4
2091  ELSE
2092  iedge = 8
2093  ENDIF
2094  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2095  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2096  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2097  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2098  ENDDO
2099 
2100  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2101  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2102  (intb==2 .OR. inte==2)) THEN
2103 
2104  dnbeg(:) = dnode(:,ijknbeg)
2105  dnend(:) = dnode(:,ijknend)
2106 
2107  arclen = 0._rfreal
2108  DO ind=indbeg+1,indend
2109  ijkn = indijk(ibeg,ind ,kbeg,inoff,ijnoff)
2110  ijkn1 = indijk(ibeg,ind-1,kbeg,inoff,ijnoff)
2111  arclen = arclen + &
2112  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2113  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2114  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2115  ENDDO
2116  ds = 0._rfreal
2117  DO ind=indbeg+1,indend-1
2118  ijkn = indijk(ibeg,ind ,kbeg,inoff,ijnoff)
2119  ijkn1 = indijk(ibeg,ind-1,kbeg,inoff,ijnoff)
2120  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2121  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2122  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2123  s = ds/arclen
2124  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2125  dnode(:,ijkn) = dn(:)
2126 ! IF (region%iRegionGlobal==10 .AND. lbound==2) &
2127 ! write(*,*)ind,s,arcLen,dNBeg,dNEnd
2128  ENDDO
2129  ENDIF ! intb
2130  ENDIF ! ipEdge
2131  ENDIF ! lbound
2132 
2133  IF (lbound==3 .OR. lbound==4) THEN
2134  IF (ipedge==1) THEN
2135  ib = ibeg
2136  ie = iend
2137  jb = jbeg
2138  je = jbeg
2139  kb = kbeg
2140  ke = kbeg
2141  indbeg = ib
2142  indend = ie
2143  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2144  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2145  intb = 1
2146  inte = 1
2147  IF (lbound==3) THEN
2148  iedge = 9
2149  ELSE
2150  iedge = 12
2151  ENDIF
2152  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2153  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2154  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2155  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2156  ENDDO
2157 
2158  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2159  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2160  (intb==2 .OR. inte==2)) THEN
2161 
2162  dnbeg(:) = dnode(:,ijknbeg)
2163  dnend(:) = dnode(:,ijknend)
2164 
2165  arclen = 0._rfreal
2166  DO ind=indbeg+1,indend
2167  ijkn = indijk(ind ,jbeg,kbeg,inoff,ijnoff)
2168  ijkn1 = indijk(ind-1,jbeg,kbeg,inoff,ijnoff)
2169  arclen = arclen + &
2170  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2171  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2172  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2173  ENDDO
2174  ds = 0._rfreal
2175  DO ind=indbeg+1,indend-1
2176  ijkn = indijk(ind ,jbeg,kbeg,inoff,ijnoff)
2177  ijkn1 = indijk(ind-1,jbeg,kbeg,inoff,ijnoff)
2178  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2179  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2180  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2181  s = ds/arclen
2182  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2183  dnode(:,ijkn) = dn(:)
2184  ENDDO
2185  ENDIF
2186  ELSEIF (ipedge==2) THEN
2187  ib = iend
2188  ie = iend
2189  jb = jbeg
2190  je = jbeg
2191  kb = kbeg
2192  ke = kend
2193  indbeg = kb
2194  indend = ke
2195  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2196  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2197  intb = 1
2198  inte = 1
2199  IF (lbound==3) THEN
2200  iedge = 5
2201  ELSE
2202  iedge = 7
2203  ENDIF
2204  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2205  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2206  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2207  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2208  ENDDO
2209 
2210  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2211  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2212  (intb==2 .OR. inte==2)) THEN
2213 
2214  dnbeg(:) = dnode(:,ijknbeg)
2215  dnend(:) = dnode(:,ijknend)
2216 
2217  arclen = 0._rfreal
2218  DO ind=indbeg+1,indend
2219  ijkn = indijk(iend,jbeg,ind ,inoff,ijnoff)
2220  ijkn1 = indijk(iend,jbeg,ind-1,inoff,ijnoff)
2221  arclen = arclen + &
2222  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2223  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2224  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2225  ENDDO
2226  ds = 0._rfreal
2227  DO ind=indbeg+1,indend-1
2228  ijkn = indijk(iend,jbeg,ind ,inoff,ijnoff)
2229  ijkn1 = indijk(iend,jbeg,ind-1,inoff,ijnoff)
2230  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2231  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2232  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2233  s = ds/arclen
2234  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2235  dnode(:,ijkn) = dn(:)
2236  ENDDO
2237  ENDIF
2238  ELSEIF (ipedge==3) THEN
2239  ib = ibeg
2240  ie = iend
2241  jb = jbeg
2242  je = jbeg
2243  kb = kend
2244  ke = kend
2245  indbeg = ib
2246  indend = ie
2247  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2248  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2249  intb = 1
2250  inte = 1
2251  IF (lbound==3) THEN
2252  iedge = 10
2253  ELSE
2254  iedge = 11
2255  ENDIF
2256  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2257  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2258  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2259  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2260  ENDDO
2261 
2262  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2263  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2264  (intb==2 .OR. inte==2)) THEN
2265 
2266  dnbeg(:) = dnode(:,ijknbeg)
2267  dnend(:) = dnode(:,ijknend)
2268 
2269  arclen = 0._rfreal
2270  DO ind=indbeg+1,indend
2271  ijkn = indijk(ind ,jbeg,kend,inoff,ijnoff)
2272  ijkn1 = indijk(ind-1,jbeg,kend,inoff,ijnoff)
2273  arclen = arclen + &
2274  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2275  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2276  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2277  ENDDO
2278  ds = 0._rfreal
2279  DO ind=indbeg+1,indend-1
2280  ijkn = indijk(ind ,jbeg,kend,inoff,ijnoff)
2281  ijkn1 = indijk(ind-1,jbeg,kend,inoff,ijnoff)
2282  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2283  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2284  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2285  s = ds/arclen
2286  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2287  dnode(:,ijkn) = dn(:)
2288  ENDDO
2289  ENDIF
2290  ELSEIF (ipedge==4) THEN
2291  ib = ibeg
2292  ie = ibeg
2293  jb = jbeg
2294  je = jbeg
2295  kb = kbeg
2296  ke = kend
2297  indbeg = kb
2298  indend = ke
2299  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2300  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2301  intb = 1
2302  inte = 1
2303  IF (lbound==3) THEN
2304  iedge = 1
2305  ELSE
2306  iedge = 3
2307  ENDIF
2308  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2309  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2310  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2311  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2312  ENDDO
2313 
2314  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2315  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2316  (intb==2 .OR. inte==2)) THEN
2317 
2318  dnbeg(:) = dnode(:,ijknbeg)
2319  dnend(:) = dnode(:,ijknend)
2320 
2321  arclen = 0._rfreal
2322  DO ind=indbeg+1,indend
2323  ijkn = indijk(ibeg,jbeg,ind ,inoff,ijnoff)
2324  ijkn1 = indijk(ibeg,jbeg,ind-1,inoff,ijnoff)
2325  arclen = arclen + &
2326  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2327  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2328  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2329  ENDDO
2330  ds = 0._rfreal
2331  DO ind=indbeg+1,indend-1
2332  ijkn = indijk(ibeg,jbeg,ind ,inoff,ijnoff)
2333  ijkn1 = indijk(ibeg,jbeg,ind-1,inoff,ijnoff)
2334  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2335  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2336  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2337  s = ds/arclen
2338  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2339  dnode(:,ijkn) = dn(:)
2340  ENDDO
2341  ENDIF ! intb
2342  ENDIF ! ipEdge
2343  ENDIF ! lbound
2344 
2345  IF (lbound==5 .OR. lbound==6) THEN
2346  IF (ipedge==1) THEN
2347  ib = ibeg
2348  ie = ibeg
2349  jb = jbeg
2350  je = jend
2351  kb = kbeg
2352  ke = kbeg
2353  indbeg = jb
2354  indend = je
2355  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2356  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2357  intb = 1
2358  inte = 1
2359  IF (lbound==5) THEN
2360  iedge = 4
2361  ELSE
2362  iedge = 2
2363  ENDIF
2364  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2365  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2366  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2367  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2368  ENDDO
2369 
2370  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2371  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2372  (intb==2 .OR. inte==2)) THEN
2373 
2374  dnbeg(:) = dnode(:,ijknbeg)
2375  dnend(:) = dnode(:,ijknend)
2376 
2377  arclen = 0._rfreal
2378  DO ind=indbeg+1,indend
2379  ijkn = indijk(ibeg,ind ,kbeg,inoff,ijnoff)
2380  ijkn1 = indijk(ibeg,ind-1,kbeg,inoff,ijnoff)
2381  arclen = arclen + &
2382  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2383  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2384  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2385  ENDDO
2386  ds = 0._rfreal
2387  DO ind=indbeg+1,indend-1
2388  ijkn = indijk(ibeg,ind ,kbeg,inoff,ijnoff)
2389  ijkn1 = indijk(ibeg,ind-1,kbeg,inoff,ijnoff)
2390  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2391  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2392  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2393  s = ds/arclen
2394  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2395  dnode(:,ijkn) = dn(:)
2396  ENDDO
2397  ENDIF
2398  ELSEIF (ipedge==2) THEN
2399  ib = ibeg
2400  ie = iend
2401  jb = jend
2402  je = jend
2403  kb = kbeg
2404  ke = kbeg
2405  indbeg = ib
2406  indend = ie
2407  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2408  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2409  intb = 1
2410  inte = 1
2411  IF (lbound==5) THEN
2412  iedge = 12
2413  ELSE
2414  iedge = 11
2415  ENDIF
2416  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2417  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2418  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2419  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2420  ENDDO
2421 
2422  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2423  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2424  (intb==2 .OR. inte==2)) THEN
2425 
2426  dnbeg(:) = dnode(:,ijknbeg)
2427  dnend(:) = dnode(:,ijknend)
2428 
2429  arclen = 0._rfreal
2430  DO ind=indbeg+1,indend
2431  ijkn = indijk(ind ,jend,kbeg,inoff,ijnoff)
2432  ijkn1 = indijk(ind-1,jend,kbeg,inoff,ijnoff)
2433  arclen = arclen + &
2434  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2435  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2436  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2437  ENDDO
2438  ds = 0._rfreal
2439  DO ind=indbeg+1,indend-1
2440  ijkn = indijk(ind ,jend,kbeg,inoff,ijnoff)
2441  ijkn1 = indijk(ind-1,jend,kbeg,inoff,ijnoff)
2442  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2443  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2444  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2445  s = ds/arclen
2446  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2447  dnode(:,ijkn) = dn(:)
2448  ENDDO
2449  ENDIF
2450  ELSEIF (ipedge==3) THEN
2451  ib = iend
2452  ie = iend
2453  jb = jbeg
2454  je = jend
2455  kb = kbeg
2456  ke = kbeg
2457  indbeg = jb
2458  indend = je
2459  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2460  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2461  intb = 1
2462  inte = 1
2463  IF (lbound==5) THEN
2464  iedge = 8
2465  ELSE
2466  iedge = 6
2467  ENDIF
2468  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2469  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2470  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2471  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2472  ENDDO
2473 
2474  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2475  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2476  (intb==2 .OR. inte==2)) THEN
2477 
2478  dnbeg(:) = dnode(:,ijknbeg)
2479  dnend(:) = dnode(:,ijknend)
2480 
2481  arclen = 0._rfreal
2482  DO ind=indbeg+1,indend
2483  ijkn = indijk(iend,ind ,kbeg,inoff,ijnoff)
2484  ijkn1 = indijk(iend,ind-1,kbeg,inoff,ijnoff)
2485  arclen = arclen + &
2486  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2487  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2488  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2489  ENDDO
2490  ds = 0._rfreal
2491  DO ind=indbeg+1,indend-1
2492  ijkn = indijk(iend,ind ,kbeg,inoff,ijnoff)
2493  ijkn1 = indijk(iend,ind-1,kbeg,inoff,ijnoff)
2494  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2495  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2496  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2497  s = ds/arclen
2498  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2499  dnode(:,ijkn) = dn(:)
2500  ENDDO
2501  ENDIF
2502  ELSEIF (ipedge==4) THEN
2503  ib = ibeg
2504  ie = iend
2505  jb = jbeg
2506  je = jbeg
2507  kb = kbeg
2508  ke = kbeg
2509  indbeg = ib
2510  indend = ie
2511  ijknbeg = indijk(ib,jb,kb,inoff,ijnoff)
2512  ijknend = indijk(ie,je,ke,inoff,ijnoff)
2513  intb = 1
2514  inte = 1
2515  IF (lbound==5) THEN
2516  iedge = 9
2517  ELSE
2518  iedge = 10
2519  ENDIF
2520  intertype = region%levels(ilev)%edgeCells(iedge)%interType
2521  DO ic = 1,grid%nCorns(region%iRegionGlobal)
2522  IF (ic==ijknbeg) intb = grid%nghbor(3,1,ic)
2523  IF (ic==ijknend) inte = grid%nghbor(3,1,ic)
2524  ENDDO
2525 
2526  IF (((intb/=0 .OR. inte/=0) .AND. (intertype==edge_interact_full .AND. &
2527  region%levels(ilev)%edgeCells(iedge)%interact .EQV. .true.)) .OR.&
2528  (intb==2 .OR. inte==2)) THEN
2529 
2530  dnbeg(:) = dnode(:,ijknbeg)
2531  dnend(:) = dnode(:,ijknend)
2532 
2533  arclen = 0._rfreal
2534  DO ind=indbeg+1,indend
2535  ijkn = indijk(ind ,jbeg,kbeg,inoff,ijnoff)
2536  ijkn1 = indijk(ind-1,jbeg,kbeg,inoff,ijnoff)
2537  arclen = arclen + &
2538  sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2539  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2540  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2541  ENDDO
2542  ds = 0._rfreal
2543  DO ind=indbeg+1,indend-1
2544  ijkn = indijk(ind ,jbeg,kbeg,inoff,ijnoff)
2545  ijkn1 = indijk(ind-1,jbeg,kbeg,inoff,ijnoff)
2546  ds = ds + sqrt((xyzold(xcoord,ijkn)-xyzold(xcoord,ijkn1))**2 + &
2547  (xyzold(ycoord,ijkn)-xyzold(ycoord,ijkn1))**2 + &
2548  (xyzold(zcoord,ijkn)-xyzold(zcoord,ijkn1))**2)
2549  s = ds/arclen
2550  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
2551  dnode(:,ijkn) = dn(:)
2552  ENDDO
2553  ENDIF ! intb
2554  ENDIF ! ipEdge
2555  ENDIF ! lbound
2556 
2557  ENDDO ! ipEdge
2558 
2559 777 CONTINUE
2560 
2561  ENDDO ! iPatch
2562 
2563 ! finalize --------------------------------------------------------------------
2564 
2565  CALL deregisterfunction( region%global )
2566 
2567 END SUBROUTINE rflo_mgframeedges
2568 
2569 !******************************************************************************
2570 !
2571 ! Purpose: exchange deformations between the regions as to ensure
2572 ! matching grid nodes at the interfaces.
2573 !
2574 ! Description: none.
2575 !
2576 ! Input: regions = data of all grid regions, deformations.
2577 !
2578 ! Output: regions%levels%grid%xyz = deformations at the boundaries.
2579 !
2580 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
2581 ! is applied to the finest grid first.
2582 !
2583 !******************************************************************************
2584 
2585 SUBROUTINE rflo_mgframeinterfaces( regions,iType )
2586 
2590  IMPLICIT NONE
2591 
2592 ! ... parameters
2593  TYPE(t_region), POINTER :: regions(:)
2594  INTEGER ::itype
2595 
2596 ! ... loop variables
2597  INTEGER :: ireg, ipatch, ipass
2598 
2599 ! ... local variables
2600  INTEGER :: bctype, iregsrc, ipatchsrc, npass
2601 
2602  TYPE(t_grid), POINTER :: grid, gridold, gridsrc
2603  TYPE(t_global), POINTER :: global
2604  TYPE(t_patch), POINTER :: patch, patchsrc
2605 
2606 !******************************************************************************
2607 
2608  global => regions(1)%global
2609 
2610  CALL registerfunction( global,'RFLO_MgFrameInterfaces',&
2611  'RFLO_ModMoveGridNconform.F90' )
2612 
2613 ! fix interfaces between regions ----------------------------------------------
2614 
2615  npass = global%moveGridNsmatch
2616  DO ipass=1,npass
2617 
2618 ! - copy / send deformations
2619 
2620  DO ireg=1,global%nRegions
2621  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
2622  regions(ireg)%active==active .AND. & ! on my processor
2623  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
2624 
2625  grid => regions(ireg)%levels(1)%grid
2626  gridold => regions(ireg)%levels(1)%gridOld
2627 
2628  DO ipatch=1,regions(ireg)%nPatches
2629  patch => regions(ireg)%levels(1)%patches(ipatch)
2630  bctype = patch%bcType
2631  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
2632  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
2633  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
2634  iregsrc = patch%srcRegion
2635  ipatchsrc = patch%srcPatch
2636  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
2637  gridsrc => regions(iregsrc)%levels(1)%grid
2638 
2639  IF (regions(iregsrc)%procid == global%myProcid) THEN
2640  CALL rflo_exchangednodecopy( regions(ireg),regions(iregsrc), &
2641  patch,patchsrc,.false., &
2642  grid%xyz,gridsrc%xyz )
2643  IF (ipass < npass) THEN
2644  CALL rflo_mgframeedgeso( regions(ireg),2,grid%boundMoved, &
2645  grid%allExternal,grid%edgeMoved, &
2646  grid%arcLen12,grid%arcLen34, &
2647  grid%arcLen56,gridold%xyzOld,grid%xyz )
2648  CALL rflo_boundarydeformation( regions(ireg), &
2649  grid%boundMoved, &
2650  grid%edgeMoved,grid%arcLen12, &
2651  grid%arcLen34,grid%arcLen56, &
2652  gridold%xyzOld,grid%xyz )
2653 
2654 ! CALL RFLO_MgFrameEdges( regions(iReg),grid%edgeMoved, &
2655 ! gridOld%xyzOld,grid%xyz )
2656 ! CALL RFLO_MgFrameBndDeformation( regions(iReg), &
2657 ! grid%boundMoved, &
2658 ! grid%edgeMoved,grid%arcLen12, &
2659 ! grid%arcLen34,grid%arcLen56, &
2660 ! gridOld%xyzOld,grid%xyz )
2661  CALL rflo_mgframerestoreexternal( regions(ireg) )
2662  ENDIF
2663  ELSE
2664  CALL rflo_exchangednodesend( regions(ireg),regions(iregsrc), &
2665  patch,grid%xyz )
2666  ENDIF
2667  ENDIF ! bcType
2668  ENDDO ! iPatch
2669 
2670  ENDIF ! region on this processor and active, grid moving
2671  ENDDO ! iReg
2672 
2673 ! - receive deformations
2674 
2675  DO ireg=1,global%nRegions
2676  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
2677  regions(ireg)%active==active .AND. & ! on my processor
2678  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
2679 
2680  grid => regions(ireg)%levels(1)%grid
2681  gridold => regions(ireg)%levels(1)%gridOld
2682 
2683  DO ipatch=1,regions(ireg)%nPatches
2684  patch => regions(ireg)%levels(1)%patches(ipatch)
2685  bctype = patch%bcType
2686  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
2687  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
2688  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
2689  iregsrc = patch%srcRegion
2690  ipatchsrc = patch%srcPatch
2691  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
2692  gridsrc => regions(iregsrc)%levels(1)%grid
2693 
2694  IF (regions(iregsrc)%procid /= global%myProcid) THEN
2695  CALL rflo_exchangednoderecv( regions(ireg),regions(iregsrc), &
2696  patch,patchsrc,.false.,grid%xyz )
2697  IF (ipass < npass) THEN
2698  CALL rflo_mgframeedgeso( regions(ireg),2,grid%boundMoved, &
2699  grid%allExternal,grid%edgeMoved, &
2700  grid%arcLen12,grid%arcLen34, &
2701  grid%arcLen56,gridold%xyzOld,grid%xyz )
2702  CALL rflo_boundarydeformation( regions(ireg), &
2703  grid%boundMoved, &
2704  grid%edgeMoved,grid%arcLen12, &
2705  grid%arcLen34,grid%arcLen56, &
2706  gridold%xyzOld,grid%xyz )
2707 
2708 ! CALL RFLO_MgFrameEdges( regions(iReg),grid%edgeMoved, &
2709 ! gridOld%xyzOld,grid%xyz )
2710 ! CALL RFLO_MgFrameBndDeformation( regions(iReg), &
2711 ! grid%boundMoved, &
2712 ! grid%edgeMoved,grid%arcLen12, &
2713 ! grid%arcLen34,grid%arcLen56, &
2714 ! gridOld%xyzOld,grid%xyz )
2715  CALL rflo_mgframerestoreexternal( regions(ireg) )
2716  ENDIF
2717 
2718  ENDIF
2719  ENDIF ! bcType
2720  ENDDO ! iPatch
2721 
2722  ENDIF ! region on this processor and active, grid moving
2723  ENDDO ! iReg
2724 
2725 ! - clear send requests
2726 
2727  DO ireg=1,global%nRegions
2728  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
2729  regions(ireg)%active==active .AND. & ! on my processor
2730  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
2731  CALL rflo_clearsendrequests( regions,ireg,.true. )
2732  ENDIF
2733  ENDDO
2734 
2735  ENDDO ! iPass
2736 
2737 ! finalize --------------------------------------------------------------------
2738 
2739  CALL deregisterfunction( global )
2740 
2741 END SUBROUTINE rflo_mgframeinterfaces
2742 
2743 !******************************************************************************
2744 !
2745 ! Purpose: calculate node displacements on those boundaries whose edges
2746 ! have moved but which were not marked as moving (finest grid only).
2747 !
2748 ! Description: none.
2749 !
2750 ! Input: region = grid dimensions
2751 ! boundMoved = flag for boundaries of a region which have moved
2752 ! edgeMoved = flag for edges whose nodes have moved
2753 ! arcLen12 = arclength between i=const. boundaries for each j, k
2754 ! arcLen34 = arclength between j=const. boundaries for each k, i
2755 ! arcLen56 = arclength between k=const. boundaries for each i, j
2756 ! xyzOld = grid from previous time step.
2757 !
2758 ! Output: dNode = updated deformations at boundaries.
2759 !
2760 ! Notes: variable dNode contains the whole 3-D field.
2761 !
2762 !******************************************************************************
2763 
2764 SUBROUTINE rflo_mgframebnddeformation( region,boundMoved,edgeMoved, &
2765  arclen12,arclen34,arclen56, &
2766  xyzold,dnode )
2767 
2769  rflo_tfint2d
2770  IMPLICIT NONE
2771 
2772 #include "Indexing.h"
2773 
2774 ! ... parameters
2775  LOGICAL :: boundmoved(6), edgemoved(12)
2776 
2777  REAL(RFREAL), POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
2778  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:)
2779 
2780  TYPE(t_region) :: region
2781 
2782 ! ... loop variables
2783  INTEGER :: ibound, l1, l2
2784 
2785 ! ... local variables
2786  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend
2787  INTEGER :: l1b, l1e, l2b, l2e, lc, ijkn, ijke(4), ijkem(4), inoff, ijnoff
2788  INTEGER :: switch(6,9)
2789 
2790  LOGICAL :: sum12
2791 
2792  REAL(RFREAL) :: arclen(4), ds(4), s(4)
2793  REAL(RFREAL) :: corner(3,8), e1(3), e2(3), e3(3), e4(3), &
2794  p1(3), p2(3), p3(3), p4(3), dn(3)
2795 
2796 !******************************************************************************
2797 
2798  CALL registerfunction( region%global,'RFLO_MgFrameBndDeformation',&
2799  'RFLO_ModMoveGridNconform.F90' )
2800 
2801 ! get dimensions --------------------------------------------------------------
2802 
2803  ilev = 1
2804  CALL rflo_getdimensphysnodes( region,ilev,ipnbeg,ipnend, &
2805  jpnbeg,jpnend,kpnbeg,kpnend )
2806  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
2807 
2808 ! set boundary switch ---------------------------------------------------------
2809 ! switch(:,1-4) = numbers of the 4 edges of a boundary
2810 ! switch(:,5-6) = first/last index in l1-direction
2811 ! switch(:,7-8) = first/last index in l2-direction
2812 ! switch(:, 9) = constant index
2813 
2814  switch(1,:) = (/ 1, 2, 3, 4, jpnbeg, jpnend, kpnbeg, kpnend, ipnbeg/)
2815  switch(2,:) = (/ 5, 6, 7, 8, jpnbeg, jpnend, kpnbeg, kpnend, ipnend/)
2816  switch(3,:) = (/ 1, 5, 9, 10, kpnbeg, kpnend, ipnbeg, ipnend, jpnbeg/)
2817  switch(4,:) = (/ 3, 7, 11, 12, kpnbeg, kpnend, ipnbeg, ipnend, jpnend/)
2818  switch(5,:) = (/ 4, 8, 9, 11, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg/)
2819  switch(6,:) = (/ 2, 6, 10, 12, ipnbeg, ipnend, jpnbeg, jpnend, kpnend/)
2820 
2821 ! store displacements at corners ----------------------------------------------
2822 
2823  corner(:,1) = dnode(:,indijk(ipnbeg,jpnbeg,kpnbeg,inoff,ijnoff))
2824  corner(:,2) = dnode(:,indijk(ipnbeg,jpnbeg,kpnend,inoff,ijnoff))
2825  corner(:,3) = dnode(:,indijk(ipnbeg,jpnend,kpnend,inoff,ijnoff))
2826  corner(:,4) = dnode(:,indijk(ipnbeg,jpnend,kpnbeg,inoff,ijnoff))
2827  corner(:,5) = dnode(:,indijk(ipnend,jpnbeg,kpnbeg,inoff,ijnoff))
2828  corner(:,6) = dnode(:,indijk(ipnend,jpnbeg,kpnend,inoff,ijnoff))
2829  corner(:,7) = dnode(:,indijk(ipnend,jpnend,kpnend,inoff,ijnoff))
2830  corner(:,8) = dnode(:,indijk(ipnend,jpnend,kpnbeg,inoff,ijnoff))
2831 
2832 ! move nodes on boundaries with active edges ----------------------------------
2833 
2834  DO ibound=1,6
2835 ! IF ((.NOT.boundMoved(iBound)) .AND. &
2836 ! (edgeMoved(switch(iBound,1)) .OR. edgeMoved(switch(iBound,2)) .OR. &
2837 ! edgeMoved(switch(iBound,3)) .OR. edgeMoved(switch(iBound,4)))) THEN
2838 
2839  IF ((edgemoved(switch(ibound,1)) .OR. edgemoved(switch(ibound,2)) .OR. &
2840  edgemoved(switch(ibound,3)) .OR. edgemoved(switch(ibound,4)))) THEN
2841 
2842  l1b = switch(ibound,5)
2843  l1e = switch(ibound,6)
2844  l2b = switch(ibound,7)
2845  l2e = switch(ibound,8)
2846  lc = switch(ibound,9)
2847 
2848  IF (ibound == 1) THEN
2849  p1(:) = corner(:,1)
2850  p2(:) = corner(:,4)
2851  p3(:) = corner(:,3)
2852  p4(:) = corner(:,2)
2853  ELSE IF (ibound == 2) THEN
2854  p1(:) = corner(:,5)
2855  p2(:) = corner(:,8)
2856  p3(:) = corner(:,7)
2857  p4(:) = corner(:,6)
2858  ELSE IF (ibound == 3) THEN
2859  p1(:) = corner(:,1)
2860  p2(:) = corner(:,2)
2861  p3(:) = corner(:,6)
2862  p4(:) = corner(:,5)
2863  ELSE IF (ibound == 4) THEN
2864  p1(:) = corner(:,4)
2865  p2(:) = corner(:,3)
2866  p3(:) = corner(:,7)
2867  p4(:) = corner(:,8)
2868  ELSE IF (ibound == 5) THEN
2869  p1(:) = corner(:,1)
2870  p2(:) = corner(:,5)
2871  p3(:) = corner(:,8)
2872  p4(:) = corner(:,4)
2873  ELSE IF (ibound == 6) THEN
2874  p1(:) = corner(:,2)
2875  p2(:) = corner(:,6)
2876  p3(:) = corner(:,7)
2877  p4(:) = corner(:,3)
2878  ENDIF
2879 
2880  ds(1:2) = 0._rfreal
2881  DO l2=l2b+1,l2e-1
2882 
2883  sum12 = .true.
2884  ds(3:4) = 0._rfreal
2885  DO l1=l1b+1,l1e-1
2886  IF (ibound==1 .OR. ibound==2) THEN
2887  ijkn = indijk(lc,l1 ,l2 ,inoff,ijnoff)
2888  ijke(1) = indijk(lc,jpnbeg,l2 ,inoff,ijnoff)
2889  ijkem(1) = indijk(lc,jpnbeg,l2-1 ,inoff,ijnoff)
2890  ijke(2) = indijk(lc,jpnend,l2 ,inoff,ijnoff)
2891  ijkem(2) = indijk(lc,jpnend,l2-1 ,inoff,ijnoff)
2892  ijke(3) = indijk(lc,l1 ,kpnbeg,inoff,ijnoff)
2893  ijkem(3) = indijk(lc,l1-1 ,kpnbeg,inoff,ijnoff)
2894  ijke(4) = indijk(lc,l1 ,kpnend,inoff,ijnoff)
2895  ijkem(4) = indijk(lc,l1-1 ,kpnend,inoff,ijnoff)
2896  arclen(1) = arclen56(lc,jpnbeg)
2897  arclen(2) = arclen56(lc,jpnend)
2898  arclen(3) = arclen34(kpnbeg,lc)
2899  arclen(4) = arclen34(kpnend,lc)
2900  ELSE IF (ibound==3 .OR. ibound==4) THEN
2901  ijkn = indijk(l2 ,lc,l1 ,inoff,ijnoff)
2902  ijke(1) = indijk(l2 ,lc,kpnbeg,inoff,ijnoff)
2903  ijkem(1) = indijk(l2-1 ,lc,kpnbeg,inoff,ijnoff)
2904  ijke(2) = indijk(l2 ,lc,kpnend,inoff,ijnoff)
2905  ijkem(2) = indijk(l2-1 ,lc,kpnend,inoff,ijnoff)
2906  ijke(3) = indijk(ipnbeg,lc,l1 ,inoff,ijnoff)
2907  ijkem(3) = indijk(ipnbeg,lc,l1-1 ,inoff,ijnoff)
2908  ijke(4) = indijk(ipnend,lc,l1 ,inoff,ijnoff)
2909  ijkem(4) = indijk(ipnend,lc,l1-1 ,inoff,ijnoff)
2910  arclen(1) = arclen12(lc,kpnbeg)
2911  arclen(2) = arclen12(lc,kpnend)
2912  arclen(3) = arclen56(ipnbeg,lc)
2913  arclen(4) = arclen56(ipnend,lc)
2914  ELSE IF (ibound==5 .OR. ibound==6) THEN
2915  ijkn = indijk(l1 ,l2 ,lc,inoff,ijnoff)
2916  ijke(1) = indijk(ipnbeg,l2 ,lc,inoff,ijnoff)
2917  ijkem(1) = indijk(ipnbeg,l2-1 ,lc,inoff,ijnoff)
2918  ijke(2) = indijk(ipnend,l2 ,lc,inoff,ijnoff)
2919  ijkem(2) = indijk(ipnend,l2-1 ,lc,inoff,ijnoff)
2920  ijke(3) = indijk(l1 ,jpnbeg,lc,inoff,ijnoff)
2921  ijkem(3) = indijk(l1-1 ,jpnbeg,lc,inoff,ijnoff)
2922  ijke(4) = indijk(l1 ,jpnend,lc,inoff,ijnoff)
2923  ijkem(4) = indijk(l1-1 ,jpnend,lc,inoff,ijnoff)
2924  arclen(1) = arclen34(lc,ipnbeg)
2925  arclen(2) = arclen34(lc,ipnend)
2926  arclen(3) = arclen12(jpnbeg,lc)
2927  arclen(4) = arclen12(jpnend,lc)
2928  ENDIF
2929  IF (sum12) THEN
2930  ds(1) = ds(1) + &
2931  sqrt((xyzold(xcoord,ijke(1))-xyzold(xcoord,ijkem(1)))**2 + &
2932  (xyzold(ycoord,ijke(1))-xyzold(ycoord,ijkem(1)))**2 + &
2933  (xyzold(zcoord,ijke(1))-xyzold(zcoord,ijkem(1)))**2)
2934  ds(2) = ds(2) + &
2935  sqrt((xyzold(xcoord,ijke(2))-xyzold(xcoord,ijkem(2)))**2 + &
2936  (xyzold(ycoord,ijke(2))-xyzold(ycoord,ijkem(2)))**2 + &
2937  (xyzold(zcoord,ijke(2))-xyzold(zcoord,ijkem(2)))**2)
2938  sum12 = .false.
2939  ENDIF
2940  ds(3) = ds(3) + &
2941  sqrt((xyzold(xcoord,ijke(3))-xyzold(xcoord,ijkem(3)))**2 + &
2942  (xyzold(ycoord,ijke(3))-xyzold(ycoord,ijkem(3)))**2 + &
2943  (xyzold(zcoord,ijke(3))-xyzold(zcoord,ijkem(3)))**2)
2944  ds(4) = ds(4) + &
2945  sqrt((xyzold(xcoord,ijke(4))-xyzold(xcoord,ijkem(4)))**2 + &
2946  (xyzold(ycoord,ijke(4))-xyzold(ycoord,ijkem(4)))**2 + &
2947  (xyzold(zcoord,ijke(4))-xyzold(zcoord,ijkem(4)))**2)
2948  s(:) = ds(:)/arclen(:)
2949  e1(:) = dnode(:,ijke(1))
2950  e2(:) = dnode(:,ijke(2))
2951  e3(:) = dnode(:,ijke(3))
2952  e4(:) = dnode(:,ijke(4))
2953  CALL rflo_tfint2d( s(1),s(2),s(3),s(4),e1,e2,e3,e4,p1,p2,p3,p4,dn )
2954  dnode(:,ijkn) = dn(:)
2955  ENDDO ! l1
2956  ENDDO ! l2
2957 
2958  ENDIF ! edgeMoved
2959  ENDDO ! iBound
2960 
2961 ! finalize --------------------------------------------------------------------
2962 
2963  CALL deregisterfunction( region%global )
2964 
2965 END SUBROUTINE rflo_mgframebnddeformation
2966 
2967 ! ******************************************************************************
2968 ! End
2969 ! ******************************************************************************
2970 
2971 END MODULE rflo_modmovegridframe
2972 
2973 ! ******************************************************************************
2974 !
2975 ! RCS Revision history:
2976 !
2977 ! $Log: RFLO_ModMoveGridNconform.F90,v $
2978 ! Revision 1.37 2009/08/27 14:04:50 mtcampbe
2979 ! Updated to enable burning motion with symmetry boundaries and enhanced
2980 ! burnout code.
2981 !
2982 ! Revision 1.36 2008/12/06 08:44:17 mtcampbe
2983 ! Updated license.
2984 !
2985 ! Revision 1.35 2008/11/19 22:17:28 mtcampbe
2986 ! Added Illinois Open Source License/Copyright
2987 !
2988 ! Revision 1.34 2006/03/18 11:02:38 wasistho
2989 ! screen printed global skewness and minvol
2990 !
2991 ! Revision 1.33 2006/03/05 22:27:41 wasistho
2992 ! fixed syntax error
2993 !
2994 ! Revision 1.32 2006/03/05 21:52:00 wasistho
2995 ! changed computational space coordinates to be based on initial grid
2996 !
2997 ! Revision 1.31 2006/02/11 03:53:18 wasistho
2998 ! made some routines public
2999 !
3000 ! Revision 1.30 2006/01/28 22:52:34 wasistho
3001 ! fixed iEdgeGlo
3002 !
3003 ! Revision 1.29 2005/11/16 22:52:27 wasistho
3004 ! update screen print moveGridFOMS
3005 !
3006 ! Revision 1.28 2005/10/28 07:41:05 wasistho
3007 ! modified FORMAT 1000
3008 !
3009 ! Revision 1.27 2005/10/27 19:19:49 wasistho
3010 ! modified screen print
3011 !
3012 ! Revision 1.26 2005/10/27 05:57:11 wasistho
3013 ! added MoveGridFoms
3014 !
3015 ! Revision 1.25 2005/08/18 19:51:12 wasistho
3016 ! added user define nPass in mgFrameInterface
3017 !
3018 ! Revision 1.24 2005/07/10 21:13:04 wasistho
3019 ! global alloc dist in SrchNeighbors and CorrectNeighbors, and added pointer grid in mgFrameMoveCorners
3020 !
3021 ! Revision 1.23 2005/06/30 19:10:27 wasistho
3022 ! made official last added conditions in mgFrameEdges
3023 !
3024 ! Revision 1.22 2005/06/29 08:38:56 wasistho
3025 ! changed RestoreExternal
3026 !
3027 ! Revision 1.21 2005/06/29 02:37:02 wasistho
3028 ! activated itype=2
3029 !
3030 ! Revision 1.20 2005/06/27 19:23:04 wasistho
3031 ! only ipass=1 in mgFrameInterfaces
3032 !
3033 ! Revision 1.19 2005/06/27 01:00:53 wasistho
3034 ! stored tolerance in tol
3035 !
3036 ! Revision 1.18 2005/06/27 00:36:52 wasistho
3037 ! change tolerance in mgFrameSrchNeighbors from 1.e-20 to 1.e-5*edgeLen
3038 !
3039 ! Revision 1.17 2005/06/26 06:25:40 wasistho
3040 ! nReg==72 to nReg==71
3041 !
3042 ! Revision 1.16 2005/06/26 06:12:00 wasistho
3043 ! adding more reegions for titan check in mgFrameSrchCorners
3044 !
3045 ! Revision 1.15 2005/06/26 05:39:44 wasistho
3046 ! bugs fixed in mgFrameCornPoints and mgFrameSrchNeighbors
3047 !
3048 ! Revision 1.14 2005/06/25 08:12:26 wasistho
3049 ! bug fixed in receiving rvar in mgFrameBCast
3050 !
3051 ! Revision 1.13 2005/06/25 06:16:15 wasistho
3052 ! swap ENDDO and ENDIF in mgframeBroadCast
3053 !
3054 ! Revision 1.12 2005/06/25 03:13:57 wasistho
3055 ! enabled nRegions /= nProcs in type 2 gridmotion
3056 !
3057 ! Revision 1.11 2005/06/19 12:36:07 wasistho
3058 ! update in mgFrameInterfaces
3059 !
3060 ! Revision 1.10 2005/06/15 19:20:40 wasistho
3061 ! simplified mgFrameInterfaces and copied RFLO_ModMoveGridNconform to RFLO_ModMoveGridFrame
3062 !
3063 ! Revision 1.9 2005/06/14 11:20:46 wasistho
3064 ! set itype to 1
3065 !
3066 ! Revision 1.8 2005/06/13 21:47:31 wasistho
3067 ! changed patch%bcCoupled to patch%bcMotion
3068 !
3069 ! Revision 1.7 2005/06/12 10:56:33 wasistho
3070 ! uncommented second/end TFI procedure
3071 !
3072 ! Revision 1.6 2005/06/12 10:12:08 wasistho
3073 ! commented second/end TFI procedure
3074 !
3075 ! Revision 1.5 2005/06/12 08:26:20 wasistho
3076 ! set initial ds=0 in mgFrameEdge
3077 !
3078 ! Revision 1.3 2005/06/12 07:59:36 wasistho
3079 ! working version of Nconform1
3080 !
3081 ! Revision 1.19 2005/06/12 06:21:29 wasistho
3082 ! modified dumax in MgFrameCorrectNeighbors
3083 !
3084 ! Revision 1.18 2005/06/12 00:50:11 wasistho
3085 ! fixed bug in defining regNc
3086 !
3087 ! Revision 1.17 2005/06/11 22:57:21 wasistho
3088 ! working version of RFLO_ModMoveGridNconform
3089 !
3090 ! Revision 1.1 2005/06/11 20:59:26 wasistho
3091 ! import as non-conforming partition version of RFLO_ModMoveGridFrame
3092 !
3093 ! Revision 1.13 2005/06/06 23:04:59 wasistho
3094 ! fixed bug put logical test between brackets
3095 !
3096 ! Revision 1.12 2005/06/05 23:03:25 wasistho
3097 ! distinguish external boundary to be fully and partly external
3098 !
3099 ! Revision 1.11 2005/06/04 22:28:06 wasistho
3100 ! more rigorous searching and made alghorithm 1 as default
3101 !
3102 ! Revision 1.10 2005/06/04 10:51:43 wasistho
3103 ! move surface containing partly external-bc
3104 !
3105 ! Revision 1.9 2005/06/04 04:25:44 wasistho
3106 ! fixed bug in determining interior vs external corners
3107 !
3108 ! Revision 1.8 2005/06/04 01:01:33 wasistho
3109 ! distinguished to AMPLIFX,Y,Z
3110 !
3111 ! Revision 1.7 2005/06/03 01:26:40 wasistho
3112 ! moved amplif to affect regCornOld gradually with iteration
3113 !
3114 ! Revision 1.6 2005/06/02 22:59:34 wasistho
3115 ! added user controlled moveGridAmplif and moveGridPower
3116 !
3117 ! Revision 1.5 2005/06/02 19:54:49 wasistho
3118 ! added control parameters amplif and pow
3119 !
3120 ! Revision 1.4 2005/06/01 22:57:33 wasistho
3121 ! commented remeshing
3122 !
3123 ! Revision 1.3 2005/06/01 08:02:39 wasistho
3124 ! increased mgFrame iteration from 5 to 10
3125 !
3126 ! Revision 1.2 2005/06/01 07:14:41 wasistho
3127 ! debuged and made more robust
3128 !
3129 ! Revision 1.1 2005/05/28 08:11:56 wasistho
3130 ! import RFLO_ModMoveGridFrame
3131 !
3132 !
3133 !
3134 ! ******************************************************************************
3135 
3136 
3137 
3138 
3139 
3140 
3141 
3142 
3143 
3144 
3145 
3146 
3147 
3148 
3149 
3150 
3151 
3152 
subroutine rflo_edgedeformation(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
**********************************************************************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 ibeg
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 rflo_calccellcentroids(region)
subroutine rflo_mgframesurfaces(regions, someMoved, iType)
j indices k indices k
Definition: Indexing.h:6
subroutine rflo_mgframeedges(region, iType, boundMoved, allExternal, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine rflo_mgframeedgeso(region, iType, boundMoved, allExternal, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
double s
Definition: blastest.C:80
**********************************************************************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 kpcbeg
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
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
NT p1
subroutine rflo_tfint2d(s1, s2, s3, s4, e1, e2, e3, e4, p1, p2, p3, p4, xyz)
Definition: RFLO_Tfint.F90:76
subroutine rflo_c2eavgcoeffs(region)
subroutine, public rflo_mgframebroadcast(regions, iselect, iter)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine rflo_extrapolategeometry(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 ic
subroutine, public rflo_movegridframe(regions)
subroutine rflo_exchangednoderecv(region, regionSrc, patch, patchSrc, average, dNode)
subroutine rflo_c2favgcoeffs(region)
double sqrt(double d)
Definition: double.h:73
**********************************************************************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 jpcbeg
subroutine rflo_changeinteriorgrid(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, 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 ipcend
subroutine rflo_movegridinterfaces(regions)
subroutine rflo_calccontrolvolumes(region)
subroutine rflo_tfint1d(s, p1, p2, xyz)
Definition: RFLO_Tfint.F90:59
**********************************************************************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_mgframemovecorners(regions)
Definition: patch.h:74
**********************************************************************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 knode iend
subroutine rflo_calcfacevectors(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 ipcbeg
subroutine rflo_exchangegeometry(regions)
subroutine rflo_movegridsurfaces(regions, someMoved)
subroutine rflo_generatecoarsegrids(region)
subroutine rflo_getpatchindicesnodes(region, patch, iLev, ibeg, iend, jbeg, jend, kbeg, kend)
logical function isnan(x)
Definition: ModTools.F90:201
subroutine rflo_mgframeinterfaces(regions, iType)
blockLoc i
Definition: read.cpp:79
**********************************************************************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
void int int REAL * x
Definition: read.cpp:74
subroutine rflo_getdeformation(region, boundMoved, dNode)
subroutine, public rflo_mgframesrchneighbors(regions)
subroutine rflo_clearsendrequests(regions, iReg, geometry)
subroutine rflo_calcgridspeeds(region)
Vector_n min(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:346
subroutine rflo_getdimensphysnodes(region, iLev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend)
j indices j
Definition: Indexing.h:6
**********************************************************************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 jpcend
**********************************************************************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 knode jend
subroutine rflo_exchangednodecopy(region, regionSrc, patch, patchSrc, average, dNode, dNodeSrc)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
**********************************************************************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 knode jbeg
subroutine rflo_boundarydeformation(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
long double dist(long double *coord1, long double *coord2, int size)
subroutine rflo_mgframerestoreexternal(region)
subroutine rflo_mgframebnddeformation(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine grid(bp)
Definition: setup_py.f90:257
**********************************************************************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 knode kbeg
subroutine deregisterfunction(global)
Definition: ModError.F90:469
**********************************************************************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_getdimensphys(region, iLev, ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend)
RT a() const
Definition: Line_2.h:140
subroutine rflo_laplacegridsmoo(regions, resid)
subroutine rflo_calcfacecentroids(region)
subroutine, public rflo_mgframecornpoints(regions)
subroutine rflo_mgframecorrectneighbors(regions)