Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLO_ModMoveGridElliptFra.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 with elliptic PDE smoothings.
26 !
27 ! Description: The procedure is similar to MOVEGRID_FRAME method, except
28 ! the Laplacian smoothing is replaced by elliptic PDE smoothing.
29 !
30 ! Notes: None.
31 !
32 ! ******************************************************************************
33 !
34 ! $Id: RFLO_ModMoveGridElliptFra.F90,v 1.16 2009/08/27 14:04:50 mtcampbe Exp $
35 !
36 ! Copyright: (c) 2004 by the University of Illinois
37 !
38 ! ******************************************************************************
39 
41 
42  USE modglobal, ONLY : t_global
43  USE moddatastruct, ONLY: t_region
44  USE modgrid, ONLY : t_grid
45  USE modbndpatch, ONLY : t_patch
46  USE modparameters
47  USE moddatatypes
48  USE moderror
49  USE modmpi
50 
51  IMPLICIT NONE
52 
53  PRIVATE
54  PUBLIC :: rflo_movegridelliptfra
55 
56 ! private : RFLO_MgElliptSurfacesFra
57 ! RFLO_MgElliptInterfacesFra
58 ! RFLO_MgElliptEdgesFra
59 ! RFLO_MgElliptBoundaries
60 
61 ! ******************************************************************************
62 ! Declarations and definitions
63 ! ******************************************************************************
64 
65  CHARACTER(CHRLEN) :: RCSIdentString = &
66  '$RCSfile: RFLO_ModMoveGridElliptFra.F90,v $ $Revision: 1.16 $'
67 
68 ! ******************************************************************************
69 ! Routines
70 ! ******************************************************************************
71 
72  CONTAINS
73 
74 
75 !******************************************************************************
76 !
77 ! Purpose: redistribute grid nodes according to the movement of the
78 ! boundaries. This function smoothes the grid globally by
79 ! volume and surface mesh smoothing based on Elliptic PDE.
80 !
81 ! Description: none.
82 !
83 ! Input: regions = data of all grid regions.
84 !
85 ! Output: regions%levels%grid%xyz = new grid coordinates.
86 !
87 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
88 ! is applied to the finest grid first.
89 !
90 !******************************************************************************
91 
92 SUBROUTINE rflo_movegridelliptfra( regions )
93 
106 
107  IMPLICIT NONE
108 
109 #ifdef GENX
110  include 'roccomf90.h'
111 #endif
112 
113 ! ... parameters
114  TYPE(t_region), POINTER :: regions(:)
115 
116 ! ... loop variables
117  INTEGER :: ireg, iter, ipatch, ijk
118 
119 ! ... local variables
120  LOGICAL :: somemoved
121 
122  INTEGER :: bctype
123 
124  REAL(RFREAL) :: resid, globalresid
125  REAL(RFREAL), POINTER :: xyz(:,:), xyzold(:,:)
126 
127  TYPE(t_grid), POINTER :: grid, gridold
128  TYPE(t_global), POINTER :: global
129  TYPE(t_patch), POINTER :: patch
130 #ifdef GENX
131  DOUBLE PRECISION :: dalpha
132 #endif
133 
134 !******************************************************************************
135 
136  global => regions(1)%global
137 
138  CALL registerfunction( global,'RFLO_MoveGridElliptFra',&
139  'RFLO_ModMoveGridElliptFra.F90' )
140 
141 #ifdef GENX
142 ! update geometry buffers -----------------------------------------------------
143 
144  dalpha = global%dtMin/global%dTimeSystem
145  CALL com_call_function( global%genxHandleGm,1,dalpha )
146 #endif
147 
148 ! receive and distribute deformations for each region -------------------------
149 
150  CALL rflo_mgelliptsurfacesfra( regions,somemoved )
151 
152 ! fix interfaces between regions ----------------------------------------------
153 
154  IF (somemoved) THEN
155  CALL rflo_mgelliptinterfacesfra( regions )
156  ENDIF
157 
158 ! update grid, dummy, corner and edge cells -----------------------------------
159 
160  DO ireg=1,global%nRegions
161  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
162  regions(ireg)%active==active .AND. & ! on my processor
163  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
164 
165 ! --- change the interior grid
166 
167  grid => regions(ireg)%levels(1)%grid
168  gridold => regions(ireg)%levels(1)%gridOld
169  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
170  grid%edgeMoved,grid%arcLen12, &
171  grid%arcLen34,grid%arcLen56, &
172  gridold%xyzOld,grid%xyz )
173 
174 ! --- update coarse grids and dummy cells
175 
176  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
177  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
178  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
179  ENDIF ! region on this processor and active, grid moving
180  ENDDO ! iReg
181  CALL rflo_exchangegeometry( regions ) ! exchange geometry
182 
183 ! smooth grid by Laplacian method ---------------------------------------------
184 
185  IF (global%moveGridNiter < 1 .AND. &
186  global%moveGridViter < 1 .AND. &
187  global%moveGridSiter < 1) THEN ! TFI only
188  IF (global%verbLevel >= verbose_high) THEN
189  IF (global%myProcid == masterproc) THEN
190  WRITE(stdout,4000) solver_name,global%skewness,global%minVol
191  WRITE(stdout,1000) solver_name, &
192  global%moveGridNiter, global%moveGridAmplifX, &
193  global%moveGridAmplifY, global%moveGridAmplifZ, &
194  global%moveGridPower, global%moveGridOrthDir, &
195  global%moveGridOrthWghtX, global%moveGridOrthWghtY, &
196  global%moveGridOrthWghtZ
197  ENDIF ! masterproc
198  ENDIF ! verbLevel
199  goto 888
200  ENDIF ! niter<1
201 
202  IF (somemoved) THEN
203  resid = 0._rfreal
204  DO iter=1,global%moveGridNiter
205  CALL rflo_laplacegridsmoo( regions,resid )
206  ENDDO
207 ! DO iter=1,global%moveGridNiter
208 ! CALL RFLO_ElliptGridSmoo( regions,resid )
209 ! ENDDO
210 
211  IF (global%verbLevel >= verbose_high) THEN
212 #ifdef MPI
213  CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
214  masterproc,global%mpiComm,global%mpierr )
215  IF (global%mpierr /= 0) CALL errorstop( global,err_mpi_trouble,&
216  __line__ )
217 #else
218  globalresid = resid
219 #endif
220  IF (global%myProcid == masterproc) THEN
221  WRITE(stdout,4000) solver_name,global%skewness,global%minVol
222  WRITE(stdout,2000) solver_name, &
223  global%moveGridNiter, global%moveGridAmplifX, &
224  global%moveGridAmplifY,global%moveGridAmplifZ, &
225  global%moveGridPower,global%moveGridOrthDir, &
226  global%moveGridOrthWghtX, global%moveGridOrthWghtY, &
227  global%moveGridOrthWghtZ, sqrt(globalresid)
228  ENDIF ! masterproc
229  ENDIF ! verbLevel
230  ENDIF ! someMoved
231 
232 ! update grid, dummy, corner and edge cells -----------------------------------
233 
234  DO ireg=1,global%nRegions
235  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
236  regions(ireg)%active==active .AND. & ! on my processor
237  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
238 
239 ! --- change xyz from coordinates to deformations
240 
241  xyz => regions(ireg)%levels(1)%grid%xyz
242  xyzold => regions(ireg)%levels(1)%gridOld%xyz
243 
244  DO ijk=lbound(xyz,2),ubound(xyz,2)
245  xyz(xcoord,ijk) = xyz(xcoord,ijk) - xyzold(xcoord,ijk)
246  xyz(ycoord,ijk) = xyz(ycoord,ijk) - xyzold(ycoord,ijk)
247  xyz(zcoord,ijk) = xyz(zcoord,ijk) - xyzold(zcoord,ijk)
248  ENDDO
249 
250 ! --- redistribute deformations at boundaries
251 
252  grid => regions(ireg)%levels(1)%grid
253  gridold => regions(ireg)%levels(1)%gridOld
254  grid%boundMoved(:) = .true.
255  grid%edgeMoved(:) = .true.
256 ! TEMPORARY - DEBUGGING OF SURFACE MOTION
257  DO ipatch=1,regions(ireg)%nPatches
258  patch => regions(ireg)%levels(1)%patches(ipatch)
259  bctype = patch%bcType
260 
261  IF (bctype .eq. bc_symmetry) THEN
262  grid%boundMoved(patch%lbound) = .false.
263  ENDIF ! bcType
264  ENDDO ! iPatch
265 ! TEMPORARY - DEBUGGING OF SURFACE MOTION
266 
267 ! WRITE(*,*) ' in here '
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  CALL rflo_mgelliptbnddeformation( regions(ireg) )
274 
275 ! --- change xyz from deformations to coordinates
276 
277  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
278  grid%edgeMoved,grid%arcLen12, &
279  grid%arcLen34,grid%arcLen56, &
280  gridold%xyzOld,grid%xyz )
281 
282 ! --- perform volume smoothing based on 3D Elliptic PDE
283 
284  IF (regions(ireg)%blockShape==region_shape_normal) THEN
285  DO iter=1,global%moveGridViter
286  CALL rflo_elliptgridsmooregion( regions(ireg),resid )
287  ENDDO
288  ENDIF
289 
290 ! --- update coarse grids and dummy cells
291 
292  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
293  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
294  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
295  ENDIF ! region on this processor and active, grid moving
296  ENDDO ! iReg
297 
298  CALL rflo_exchangegeometry( regions ) ! exchange geometry
299 
300 ! print residual of 3D Elliptic PDE smoothing ---------------------------------
301 
302  IF (global%verbLevel >= verbose_high) THEN
303 #ifdef MPI
304  CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
305  masterproc,global%mpiComm,global%mpierr )
306  IF (global%mpierr /= 0) CALL errorstop( global,err_mpi_trouble,&
307  __line__ )
308 #else
309  globalresid = resid
310 #endif
311  IF (global%myProcid == masterproc) THEN
312  WRITE(stdout,3000) solver_name,global%moveGridViter,sqrt(globalresid)
313  ENDIF
314  ENDIF ! verbLevel
315 
316 888 CONTINUE
317 
318 ! calculate new metrics and grid speeds ---------------------------------------
319 
320  DO ireg=1,global%nRegions
321  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
322  regions(ireg)%active==active .AND. & ! on my processor
323  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
324  CALL rflo_calcfacevectors( regions(ireg) ) ! faces
325  CALL rflo_calccontrolvolumes( regions(ireg) ) ! volumes
326  CALL rflo_calccellcentroids( regions(ireg) ) ! cell centroids
327  IF (regions(ireg)%mixtInput%faceEdgeAvg==fe_avg_linear) &
328  CALL rflo_c2favgcoeffs( regions(ireg) ) ! cell2face averaging
329  CALL rflo_c2eavgcoeffs( regions(ireg) ) ! cell2edge averaging
330  CALL rflo_checkmetrics( ireg,regions(ireg) ) ! check metrics
331  CALL rflo_calcgridspeeds( regions(ireg) ) ! grid speeds
332  ENDIF ! region on this processor and active, grid moving
333  ENDDO ! iReg
334 
335 ! global grid quality measure -------------------------------------------------
336 
337  CALL rflo_gridqualityglobal( regions )
338 
339 ! finalize --------------------------------------------------------------------
340 
341  CALL deregisterfunction( global )
342 
343 1000 FORMAT(a,1x,'Global-Elliptic-PDE grid motion:', &
344  i5,1x,4(1pe9.2),i4,3(1pe9.2))
345 2000 FORMAT(a,1x,'Global-Elliptic-PDE grid motion:', &
346  i5,1x,4(1pe9.2),i4,3(1pe9.2),1pe12.4)
347 3000 FORMAT(a,1x,'Global-Elliptic-PDE grid motion:', &
348  i5,1pe12.4)
349 4000 FORMAT(a,1x,'global skewness, minvol:',2(1pe14.5))
350 
351 END SUBROUTINE rflo_movegridelliptfra
352 
353 
354 !******************************************************************************
355 !
356 ! Purpose: receive and distribute the deformations of surfaces
357 ! in block-wise manner.
358 !
359 ! Description: surface smoothing is based on TFI
360 !
361 ! Input: regions = data of all grid regions.
362 !
363 ! Output: regions%levels%grid%xyz = deformations at the boundaries
364 ! someMoved = parts of grid moved.
365 !
366 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
367 ! is applied to the finest grid first.
368 !
369 !******************************************************************************
370 
371 SUBROUTINE rflo_mgelliptsurfacesfra( regions,someMoved )
372 
380  IMPLICIT NONE
381 
382 ! ... parameters
383  LOGICAL :: somemoved
384 
385  TYPE(t_region), POINTER :: regions(:)
386 
387 ! ... loop variables
388  INTEGER :: ireg, iter
389 
390 ! ... local variables
391  TYPE(t_grid), POINTER :: grid, gridold
392  TYPE(t_global), POINTER :: global
393 
394 !******************************************************************************
395 
396  global => regions(1)%global
397 
398  CALL registerfunction( global,'RFLO_MgElliptSurfacesFra',&
399  'RFLO_ModMoveGridElliptFra.F90' )
400 
401 ! obtain external deformations ------------------------------------------------
402 
403  somemoved = .false.
404 
405  DO ireg=1,global%nRegions
406  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
407  regions(ireg)%active==active .AND. & ! on my processor
408  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
409 
410  grid => regions(ireg)%levels(1)%grid
411  gridold => regions(ireg)%levels(1)%gridOld
412  somemoved = .true.
413 
414 ! --- store the old grid
415 
416  gridold%indSvel = grid%indSvel
417  gridold%ipc = grid%ipc
418  gridold%jpc = grid%jpc
419  gridold%kpc = grid%kpc
420  gridold%xyz(:,:) = grid%xyz(:,:)
421  gridold%si(:,:) = grid%si(:,:)
422  gridold%sj(:,:) = grid%sj(:,:)
423  gridold%sk(:,:) = grid%sk(:,:)
424  gridold%vol(:) = grid%vol(:)
425 
426 ! --- calculate arclengths between boundaries
427 
428  CALL rflo_arclengthbounds( regions(ireg),gridold%xyzOld, &
429  grid%arcLen12,grid%arcLen34,grid%arcLen56 )
430 
431 ! --- get the boundary deformations
432 
433  CALL rflo_getdeformation( regions(ireg),grid%boundMoved,grid%xyz )
434 
435  grid%xyzOld(:,:) = grid%xyz(:,:)
436 
437  ENDIF ! region on this processor and active, grid moving
438  ENDDO ! iReg
439 
440 ! broadcast and compute block corners deformation -----------------------------
441 
442  iter = 1
443  CALL rflo_mgframebroadcast( regions,1,iter )
444  CALL rflo_mgframecorrectneighbors( regions )
445  CALL rflo_mgframemovecorners( regions )
446 
447  DO iter = 2,10
448  CALL rflo_mgframebroadcast( regions,1,iter )
449  CALL rflo_mgframemovecorners( regions )
450  ENDDO
451  CALL rflo_mgframebroadcast( regions,2,1 ) ! broadcast cBuff
452  CALL rflo_mgframeorthoshift( regions ) ! orth. to solid surf.
453 
454 ! compute edge and boundary deformations --------------------------------------
455 
456  DO ireg=1,global%nRegions
457  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
458  regions(ireg)%active==active .AND. & ! on my processor
459  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
460 
461  grid => regions(ireg)%levels(1)%grid
462  gridold => regions(ireg)%levels(1)%gridOld
463 
464 ! --- calculate deformations at remaining edges
465 
466  grid%edgeMoved(:) = .false.
467  CALL rflo_mgelliptedgesfra( regions(ireg),grid%boundMoved, &
468  grid%allExternal,grid%edgeMoved,grid%arcLen12, &
469  grid%arcLen34,grid%arcLen56,gridold%xyzOld, &
470  gridold%xyz,grid%xyz )
471 
472 ! --- correct deformations at straight edges
473 
474  IF (global%moveGridNiter < 1) THEN
475  CALL rflo_edgedeformationstraight( regions(ireg),grid%boundMoved, &
476  grid%edgeStraight,grid%edgeMoved, &
477  grid%arcLen12,grid%arcLen34,grid%arcLen56, &
478  gridold%xyzOld,gridold%xyz,grid%xyz )
479  ENDIF
480 
481 ! --- calculate deformations at remaining boundaries
482 
483  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
484  grid%edgeMoved,grid%arcLen12, &
485  grid%arcLen34,grid%arcLen56, &
486  gridold%xyzOld,grid%xyz )
487 
488  ENDIF ! region on this processor and active, grid moving
489  ENDDO ! iReg
490 
491 ! finalize --------------------------------------------------------------------
492 
493  CALL deregisterfunction( global )
494 
495 END SUBROUTINE rflo_mgelliptsurfacesfra
496 
497 
498 !******************************************************************************
499 !
500 ! Purpose: calculate node displacements on those edges whose end points have
501 ! moved, but the associated boundaries were not updated yet (finest
502 ! grid only).
503 !
504 ! Description: points along an edge are shifted using 1-D linear transfinite
505 ! interpolation (TFI).
506 !
507 ! Input: region = grid dimensions
508 ! boundMoved = flag for boundaries of a region which have moved
509 ! arcLen12 = arclength between i=const. boundaries for each j, k
510 ! arcLen34 = arclength between j=const. boundaries for each k, i
511 ! arcLen56 = arclength between k=const. boundaries for each i, j
512 ! xyzOld = initial grid.
513 !
514 ! Output: edgeMoved = flag if discretization at an edge was changed
515 ! dNode = updated deformations at edges.
516 !
517 ! Notes: variable dNode contains the whole 3-D field.
518 !
519 !******************************************************************************
520 
521 SUBROUTINE rflo_mgelliptedgesfra( region,boundMoved,allExternal,edgeMoved, &
522  arclen12,arclen34,arclen56,xyzref,xyzold, &
523  dnode )
524 
527  IMPLICIT NONE
528 
529 #include "Indexing.h"
530 
531 ! ... parameters
532  LOGICAL :: boundmoved(6), allexternal(6), edgemoved(12)
533 
534  REAL(RFREAL), POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
535  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:), xyzref(:,:)
536 
537  TYPE(t_region) :: region
538 
539 ! ... loop variables
540  INTEGER :: iedge, ind
541 
542 ! ... local variables
543  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend, l1c, l2c
544  INTEGER :: indbeg, indend, ijkn, ijkn1, ijknbeg, ijknend, inoff, ijnoff
545  INTEGER :: switch(12,11), intertype, iedgeglo
546 
547  REAL(RFREAL) :: arclen, ds, s, dn(3), dnbeg(3), dnend(3)
548  LOGICAL :: interact
549 
550 !******************************************************************************
551 
552  CALL registerfunction( region%global,'RFLO_MgElliptEdgesFra',&
553  'RFLO_ModMoveGridElliptFra.F90')
554 
555 ! get dimensions --------------------------------------------------------------
556 
557  ilev = 1
558  CALL rflo_getdimensphysnodes( region,ilev,ipnbeg,ipnend, &
559  jpnbeg,jpnend,kpnbeg,kpnend )
560  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
561 
562 ! set edge switch -------------------------------------------------------------
563 ! switch(:,1) = begins at boundary
564 ! switch(:,2) = ends on boundary
565 ! switch(:,3) = right boundary
566 ! switch(:,4) = left boundary
567 ! switch(:,5) = direction (from-to boundary)
568 ! switch(:,6) = start index
569 ! switch(:,7) = end index
570 ! switch(:,8) = constant index in 1st direction
571 ! switch(:,9) = constant index in 2nd direction
572 ! switch(:,10) = start corner number
573 ! switch(:,11) = end corner number
574 
575  switch( 1,:) = (/5, 6, 1, 3, 56, kpnbeg, kpnend, ipnbeg, jpnbeg, 1, 2/)
576  switch( 2,:) = (/3, 4, 1, 6, 34, jpnbeg, jpnend, kpnend, ipnbeg, 2, 3/)
577  switch( 3,:) = (/5, 6, 1, 4, 56, kpnbeg, kpnend, ipnbeg, jpnend, 4, 3/)
578  switch( 4,:) = (/3, 4, 1, 5, 34, jpnbeg, jpnend, kpnbeg, ipnbeg, 1, 4/)
579  switch( 5,:) = (/5, 6, 2, 3, 56, kpnbeg, kpnend, ipnend, jpnbeg, 5, 6/)
580  switch( 6,:) = (/3, 4, 2, 6, 34, jpnbeg, jpnend, kpnend, ipnend, 6, 7/)
581  switch( 7,:) = (/5, 6, 2, 4, 56, kpnbeg, kpnend, ipnend, jpnend, 8, 7/)
582  switch( 8,:) = (/3, 4, 2, 5, 34, jpnbeg, jpnend, kpnbeg, ipnend, 5, 8/)
583  switch( 9,:) = (/1, 2, 3, 5, 12, ipnbeg, ipnend, jpnbeg, kpnbeg, 1, 5/)
584  switch(10,:) = (/1, 2, 3, 6, 12, ipnbeg, ipnend, jpnbeg, kpnend, 2, 6/)
585  switch(11,:) = (/1, 2, 4, 5, 12, ipnbeg, ipnend, jpnend, kpnbeg, 4, 8/)
586  switch(12,:) = (/1, 2, 4, 6, 12, ipnbeg, ipnend, jpnend, kpnend, 3, 7/)
587 
588 ! edge movement flag ----------------------------------------------------------
589 
590  IF (boundmoved(1)) THEN
591  edgemoved( 1) = .true.; edgemoved( 2) = .true.
592  edgemoved( 3) = .true.; edgemoved( 4) = .true.
593  ENDIF
594  IF (boundmoved(2)) THEN
595  edgemoved( 5) = .true.; edgemoved( 6) = .true.
596  edgemoved( 7) = .true.; edgemoved( 8) = .true.
597  ENDIF
598  IF (boundmoved(3)) THEN
599  edgemoved( 1) = .true.; edgemoved( 5) = .true.
600  edgemoved( 9) = .true.; edgemoved(10) = .true.
601  ENDIF
602  IF (boundmoved(4)) THEN
603  edgemoved( 3) = .true.; edgemoved( 7) = .true.
604  edgemoved(11) = .true.; edgemoved(12) = .true.
605  ENDIF
606  IF (boundmoved(5)) THEN
607  edgemoved( 4) = .true.; edgemoved( 8) = .true.
608  edgemoved( 9) = .true.; edgemoved(11) = .true.
609  ENDIF
610  IF (boundmoved(6)) THEN
611  edgemoved( 2) = .true.; edgemoved( 6) = .true.
612  edgemoved(10) = .true.; edgemoved(12) = .true.
613  ENDIF
614 
615 ! loop over all 12 edges ------------------------------------------------------
616 
617  DO iedge=1,12
618  IF (.NOT.edgemoved(iedge)) THEN
619 
620  edgemoved(iedge) = .true.
621 
622  ds = 0._rfreal
623  indbeg = switch(iedge,6)
624  indend = switch(iedge,7)
625  l1c = switch(iedge,8)
626  l2c = switch(iedge,9)
627 
628  iedgeglo = iedge
629  IF (iedge==11) iedgeglo=12
630  IF (iedge==12) iedgeglo=11
631  interact = region%levels(ilev)%edgeCells(iedgeglo)%interact
632  intertype = region%levels(ilev)%edgeCells(iedgeglo)%interType
633 
634  IF (((region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==1 .OR. &
635  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==1) .AND. &
636  ((interact .EQV. .true.) .AND. (intertype==edge_interact_full))) &
637  .OR. &
638  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,10))==2 .OR. &
639  region%levels(ilev)%grid%nghbor(3,1,switch(iedge,11))==2) THEN
640 
641 ! IF (region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,10))==1 .OR. &
642 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,11))==1 .OR. &
643 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,10))==2 .OR. &
644 ! region%levels(iLev)%grid%nghbor(3,1,switch(iEdge,11))==2) THEN
645 
646  DO ind=indbeg+1,indend-1
647  IF (switch(iedge,5) == 12) THEN
648  ijkn = indijk(ind ,l1c,l2c,inoff,ijnoff)
649  ijkn1 = indijk(ind-1 ,l1c,l2c,inoff,ijnoff)
650  ijknbeg = indijk(indbeg,l1c,l2c,inoff,ijnoff)
651  ijknend = indijk(indend,l1c,l2c,inoff,ijnoff)
652  arclen = arclen12(l1c,l2c)
653  ELSE IF (switch(iedge,5) == 34) THEN
654  ijkn = indijk(l2c,ind ,l1c,inoff,ijnoff)
655  ijkn1 = indijk(l2c,ind-1 ,l1c,inoff,ijnoff)
656  ijknbeg = indijk(l2c,indbeg,l1c,inoff,ijnoff)
657  ijknend = indijk(l2c,indend,l1c,inoff,ijnoff)
658  arclen = arclen34(l1c,l2c)
659  ELSE IF (switch(iedge,5) == 56) THEN
660  ijkn = indijk(l1c,l2c,ind ,inoff,ijnoff)
661  ijkn1 = indijk(l1c,l2c,ind-1 ,inoff,ijnoff)
662  ijknbeg = indijk(l1c,l2c,indbeg,inoff,ijnoff)
663  ijknend = indijk(l1c,l2c,indend,inoff,ijnoff)
664  arclen = arclen56(l1c,l2c)
665  ENDIF
666  dnbeg(:) = dnode(:,ijknbeg) + xyzold(:,ijknbeg) - xyzref(:,ijknbeg)
667  dnend(:) = dnode(:,ijknend) + xyzold(:,ijknend) - xyzref(:,ijknend)
668 
669  ds = ds + sqrt((xyzref(xcoord,ijkn)-xyzref(xcoord,ijkn1))**2 + &
670  (xyzref(ycoord,ijkn)-xyzref(ycoord,ijkn1))**2 + &
671  (xyzref(zcoord,ijkn)-xyzref(zcoord,ijkn1))**2)
672  s = ds/arclen
673 
674  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
675  dnode(:,ijkn) = dn(:) + xyzref(:,ijkn) - xyzold(:,ijkn)
676  ENDDO ! i
677  ENDIF ! nghbor
678  ENDIF ! edgeMoved
679  ENDDO ! iEdge
680 
681 ! finalize --------------------------------------------------------------------
682 
683  CALL deregisterfunction( region%global )
684 
685 END SUBROUTINE rflo_mgelliptedgesfra
686 
687 
688 !******************************************************************************
689 !
690 ! Purpose: calculate node displacements on those boundaries whose edges
691 ! have moved but which were not marked as moving (finest grid only).
692 !
693 ! Description: none.
694 !
695 ! Input: region = grid dimensions
696 ! boundMoved = flag for boundaries of a region which have moved
697 ! edgeMoved = flag for edges whose nodes have moved
698 ! arcLen12 = arclength between i=const. boundaries for each j, k
699 ! arcLen34 = arclength between j=const. boundaries for each k, i
700 ! arcLen56 = arclength between k=const. boundaries for each i, j
701 ! xyzOld = initial grid.
702 !
703 ! Output: dNode = updated deformations at boundaries.
704 !
705 ! Notes: variable dNode contains the whole 3-D field.
706 !
707 !******************************************************************************
708 
709 SUBROUTINE rflo_mgelliptboundaries( region,boundMoved,edgeMoved, &
710  arclen12,arclen34,arclen56, &
711  xyzold,dnode )
712 
715  IMPLICIT NONE
716 
717 #include "Indexing.h"
718 
719 ! ... parameters
720  LOGICAL :: boundmoved(6), edgemoved(12)
721 
722  REAL(RFREAL), POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
723  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:)
724 
725  TYPE(t_region) :: region
726 
727 ! ... loop variables
728  INTEGER :: ibound, l1, l2
729 
730 ! ... local variables
731  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend
732  INTEGER :: l1b, l1e, l2b, l2e, lc, ijkn, ijke(4), ijkem(4), inoff, ijnoff
733  INTEGER :: switch(6,9)
734 
735  LOGICAL :: sum12
736 
737  REAL(RFREAL) :: arclen(4), ds(4), s(4)
738  REAL(RFREAL) :: corner(3,8), e1(3), e2(3), e3(3), e4(3), &
739  p1(3), p2(3), p3(3), p4(3), dn(3)
740 
741 !******************************************************************************
742 
743  CALL registerfunction( region%global,'RFLO_MgElliptBoundaries',&
744  'RFLO_ModMoveGridElliptFra.F90' )
745 
746 ! get dimensions --------------------------------------------------------------
747 
748  ilev = 1
749  CALL rflo_getdimensphysnodes( region,ilev,ipnbeg,ipnend, &
750  jpnbeg,jpnend,kpnbeg,kpnend )
751  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
752 
753 ! set boundary switch ---------------------------------------------------------
754 ! switch(:,1-4) = numbers of the 4 edges of a boundary
755 ! switch(:,5-6) = first/last index in l1-direction
756 ! switch(:,7-8) = first/last index in l2-direction
757 ! switch(:, 9) = constant index
758 
759  switch(1,:) = (/ 1, 2, 3, 4, jpnbeg, jpnend, kpnbeg, kpnend, ipnbeg/)
760  switch(2,:) = (/ 5, 6, 7, 8, jpnbeg, jpnend, kpnbeg, kpnend, ipnend/)
761  switch(3,:) = (/ 1, 5, 9, 10, kpnbeg, kpnend, ipnbeg, ipnend, jpnbeg/)
762  switch(4,:) = (/ 3, 7, 11, 12, kpnbeg, kpnend, ipnbeg, ipnend, jpnend/)
763  switch(5,:) = (/ 4, 8, 9, 11, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg/)
764  switch(6,:) = (/ 2, 6, 10, 12, ipnbeg, ipnend, jpnbeg, jpnend, kpnend/)
765 
766 ! store displacements at corners ----------------------------------------------
767 
768  corner(:,1) = dnode(:,indijk(ipnbeg,jpnbeg,kpnbeg,inoff,ijnoff))
769  corner(:,2) = dnode(:,indijk(ipnbeg,jpnbeg,kpnend,inoff,ijnoff))
770  corner(:,3) = dnode(:,indijk(ipnbeg,jpnend,kpnend,inoff,ijnoff))
771  corner(:,4) = dnode(:,indijk(ipnbeg,jpnend,kpnbeg,inoff,ijnoff))
772  corner(:,5) = dnode(:,indijk(ipnend,jpnbeg,kpnbeg,inoff,ijnoff))
773  corner(:,6) = dnode(:,indijk(ipnend,jpnbeg,kpnend,inoff,ijnoff))
774  corner(:,7) = dnode(:,indijk(ipnend,jpnend,kpnend,inoff,ijnoff))
775  corner(:,8) = dnode(:,indijk(ipnend,jpnend,kpnbeg,inoff,ijnoff))
776 
777 ! move nodes on boundaries with active edges ----------------------------------
778 
779  DO ibound=1,6
780  IF ((.NOT.boundmoved(ibound)) .AND. &
781  (edgemoved(switch(ibound,1)) .OR. edgemoved(switch(ibound,2)) .OR. &
782  edgemoved(switch(ibound,3)) .OR. edgemoved(switch(ibound,4)))) THEN
783 
784 ! IF ((edgeMoved(switch(iBound,1)) .OR. edgeMoved(switch(iBound,2)) .OR. &
785 ! edgeMoved(switch(iBound,3)) .OR. edgeMoved(switch(iBound,4)))) THEN
786 
787  l1b = switch(ibound,5)
788  l1e = switch(ibound,6)
789  l2b = switch(ibound,7)
790  l2e = switch(ibound,8)
791  lc = switch(ibound,9)
792 
793  IF (ibound == 1) THEN
794  p1(:) = corner(:,1)
795  p2(:) = corner(:,4)
796  p3(:) = corner(:,3)
797  p4(:) = corner(:,2)
798  ELSE IF (ibound == 2) THEN
799  p1(:) = corner(:,5)
800  p2(:) = corner(:,8)
801  p3(:) = corner(:,7)
802  p4(:) = corner(:,6)
803  ELSE IF (ibound == 3) THEN
804  p1(:) = corner(:,1)
805  p2(:) = corner(:,2)
806  p3(:) = corner(:,6)
807  p4(:) = corner(:,5)
808  ELSE IF (ibound == 4) THEN
809  p1(:) = corner(:,4)
810  p2(:) = corner(:,3)
811  p3(:) = corner(:,7)
812  p4(:) = corner(:,8)
813  ELSE IF (ibound == 5) THEN
814  p1(:) = corner(:,1)
815  p2(:) = corner(:,5)
816  p3(:) = corner(:,8)
817  p4(:) = corner(:,4)
818  ELSE IF (ibound == 6) THEN
819  p1(:) = corner(:,2)
820  p2(:) = corner(:,6)
821  p3(:) = corner(:,7)
822  p4(:) = corner(:,3)
823  ENDIF
824 
825  ds(1:2) = 0._rfreal
826  DO l2=l2b+1,l2e-1
827 
828  sum12 = .true.
829  ds(3:4) = 0._rfreal
830  DO l1=l1b+1,l1e-1
831  IF (ibound==1 .OR. ibound==2) THEN
832  ijkn = indijk(lc,l1 ,l2 ,inoff,ijnoff)
833  ijke(1) = indijk(lc,jpnbeg,l2 ,inoff,ijnoff)
834  ijkem(1) = indijk(lc,jpnbeg,l2-1 ,inoff,ijnoff)
835  ijke(2) = indijk(lc,jpnend,l2 ,inoff,ijnoff)
836  ijkem(2) = indijk(lc,jpnend,l2-1 ,inoff,ijnoff)
837  ijke(3) = indijk(lc,l1 ,kpnbeg,inoff,ijnoff)
838  ijkem(3) = indijk(lc,l1-1 ,kpnbeg,inoff,ijnoff)
839  ijke(4) = indijk(lc,l1 ,kpnend,inoff,ijnoff)
840  ijkem(4) = indijk(lc,l1-1 ,kpnend,inoff,ijnoff)
841  arclen(1) = arclen56(lc,jpnbeg)
842  arclen(2) = arclen56(lc,jpnend)
843  arclen(3) = arclen34(kpnbeg,lc)
844  arclen(4) = arclen34(kpnend,lc)
845  ELSE IF (ibound==3 .OR. ibound==4) THEN
846  ijkn = indijk(l2 ,lc,l1 ,inoff,ijnoff)
847  ijke(1) = indijk(l2 ,lc,kpnbeg,inoff,ijnoff)
848  ijkem(1) = indijk(l2-1 ,lc,kpnbeg,inoff,ijnoff)
849  ijke(2) = indijk(l2 ,lc,kpnend,inoff,ijnoff)
850  ijkem(2) = indijk(l2-1 ,lc,kpnend,inoff,ijnoff)
851  ijke(3) = indijk(ipnbeg,lc,l1 ,inoff,ijnoff)
852  ijkem(3) = indijk(ipnbeg,lc,l1-1 ,inoff,ijnoff)
853  ijke(4) = indijk(ipnend,lc,l1 ,inoff,ijnoff)
854  ijkem(4) = indijk(ipnend,lc,l1-1 ,inoff,ijnoff)
855  arclen(1) = arclen12(lc,kpnbeg)
856  arclen(2) = arclen12(lc,kpnend)
857  arclen(3) = arclen56(ipnbeg,lc)
858  arclen(4) = arclen56(ipnend,lc)
859  ELSE IF (ibound==5 .OR. ibound==6) THEN
860  ijkn = indijk(l1 ,l2 ,lc,inoff,ijnoff)
861  ijke(1) = indijk(ipnbeg,l2 ,lc,inoff,ijnoff)
862  ijkem(1) = indijk(ipnbeg,l2-1 ,lc,inoff,ijnoff)
863  ijke(2) = indijk(ipnend,l2 ,lc,inoff,ijnoff)
864  ijkem(2) = indijk(ipnend,l2-1 ,lc,inoff,ijnoff)
865  ijke(3) = indijk(l1 ,jpnbeg,lc,inoff,ijnoff)
866  ijkem(3) = indijk(l1-1 ,jpnbeg,lc,inoff,ijnoff)
867  ijke(4) = indijk(l1 ,jpnend,lc,inoff,ijnoff)
868  ijkem(4) = indijk(l1-1 ,jpnend,lc,inoff,ijnoff)
869  arclen(1) = arclen34(lc,ipnbeg)
870  arclen(2) = arclen34(lc,ipnend)
871  arclen(3) = arclen12(jpnbeg,lc)
872  arclen(4) = arclen12(jpnend,lc)
873  ENDIF
874  IF (sum12) THEN
875  ds(1) = ds(1) + &
876  sqrt((xyzold(xcoord,ijke(1))-xyzold(xcoord,ijkem(1)))**2 + &
877  (xyzold(ycoord,ijke(1))-xyzold(ycoord,ijkem(1)))**2 + &
878  (xyzold(zcoord,ijke(1))-xyzold(zcoord,ijkem(1)))**2)
879  ds(2) = ds(2) + &
880  sqrt((xyzold(xcoord,ijke(2))-xyzold(xcoord,ijkem(2)))**2 + &
881  (xyzold(ycoord,ijke(2))-xyzold(ycoord,ijkem(2)))**2 + &
882  (xyzold(zcoord,ijke(2))-xyzold(zcoord,ijkem(2)))**2)
883  sum12 = .false.
884  ENDIF
885  ds(3) = ds(3) + &
886  sqrt((xyzold(xcoord,ijke(3))-xyzold(xcoord,ijkem(3)))**2 + &
887  (xyzold(ycoord,ijke(3))-xyzold(ycoord,ijkem(3)))**2 + &
888  (xyzold(zcoord,ijke(3))-xyzold(zcoord,ijkem(3)))**2)
889  ds(4) = ds(4) + &
890  sqrt((xyzold(xcoord,ijke(4))-xyzold(xcoord,ijkem(4)))**2 + &
891  (xyzold(ycoord,ijke(4))-xyzold(ycoord,ijkem(4)))**2 + &
892  (xyzold(zcoord,ijke(4))-xyzold(zcoord,ijkem(4)))**2)
893  s(:) = ds(:)/arclen(:)
894  e1(:) = dnode(:,ijke(1))
895  e2(:) = dnode(:,ijke(2))
896  e3(:) = dnode(:,ijke(3))
897  e4(:) = dnode(:,ijke(4))
898  CALL rflo_tfint2d( s(1),s(2),s(3),s(4),e1,e2,e3,e4,p1,p2,p3,p4,dn )
899  dnode(:,ijkn) = dn(:)
900  ENDDO ! l1
901  ENDDO ! l2
902 
903  ENDIF ! edgeMoved
904  ENDDO ! iBound
905 
906 ! finalize --------------------------------------------------------------------
907 
908  CALL deregisterfunction( region%global )
909 
910 END SUBROUTINE rflo_mgelliptboundaries
911 
912 
913 !******************************************************************************
914 !
915 ! Purpose: exchange deformations between the regions as to ensure
916 ! matching grid nodes at the interfaces.
917 !
918 ! Description: none.
919 !
920 ! Input: regions = data of all grid regions, deformations.
921 !
922 ! Output: regions%levels%grid%xyz = deformations at the boundaries.
923 !
924 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
925 ! is applied to the finest grid first.
926 !
927 !******************************************************************************
928 
929 SUBROUTINE rflo_mgelliptinterfacesfra( regions )
930 
934  IMPLICIT NONE
935 
936 ! ... parameters
937  TYPE(t_region), POINTER :: regions(:)
938 
939 ! ... loop variables
940  INTEGER :: ireg, ipatch, ipass
941 
942 ! ... local variables
943  INTEGER :: bctype, iregsrc, ipatchsrc
944 
945  TYPE(t_grid), POINTER :: grid, gridold, gridsrc
946  TYPE(t_global), POINTER :: global
947  TYPE(t_patch), POINTER :: patch, patchsrc
948 
949 !******************************************************************************
950 
951  global => regions(1)%global
952 
953  CALL registerfunction( global,'RFLO_MgElliptInterfacesFra',&
954  'RFLO_ModMoveGridElliptFra.F90' )
955 
956 ! fix interfaces between regions ----------------------------------------------
957 
958  DO ipass=1,2
959 
960 ! - copy / send deformations
961 
962  DO ireg=1,global%nRegions
963  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
964  regions(ireg)%active==active .AND. & ! on my processor
965  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
966 
967  grid => regions(ireg)%levels(1)%grid
968  gridold => regions(ireg)%levels(1)%gridOld
969 
970  DO ipatch=1,regions(ireg)%nPatches
971  patch => regions(ireg)%levels(1)%patches(ipatch)
972  bctype = patch%bcType
973  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
974  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
975  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
976  iregsrc = patch%srcRegion
977  ipatchsrc = patch%srcPatch
978  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
979  gridsrc => regions(iregsrc)%levels(1)%grid
980 
981  IF (regions(iregsrc)%procid == global%myProcid) THEN
982  CALL rflo_exchangednodecopy( regions(ireg),regions(iregsrc), &
983  patch,patchsrc,.false., &
984  grid%xyz,gridsrc%xyz )
985  CALL rflo_mgelliptedgesfra( regions(ireg),grid%boundMoved, &
986  grid%allExternal,grid%edgeMoved, &
987  grid%arcLen12,grid%arcLen34, &
988  grid%arcLen56,gridold%xyzOld, &
989  gridold%xyz,grid%xyz )
990  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
991  grid%edgeMoved,grid%arcLen12, &
992  grid%arcLen34,grid%arcLen56, &
993  gridold%xyzOld,grid%xyz )
994  ELSE
995  CALL rflo_exchangednodesend( regions(ireg),regions(iregsrc), &
996  patch,grid%xyz )
997  ENDIF
998  ENDIF ! bcType
999  ENDDO ! iPatch
1000 
1001  ENDIF ! region on this processor and active, grid moving
1002  ENDDO ! iReg
1003 
1004 ! - receive deformations
1005 
1006  DO ireg=1,global%nRegions
1007  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1008  regions(ireg)%active==active .AND. & ! on my processor
1009  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
1010 
1011  grid => regions(ireg)%levels(1)%grid
1012  gridold => regions(ireg)%levels(1)%gridOld
1013 
1014  DO ipatch=1,regions(ireg)%nPatches
1015  patch => regions(ireg)%levels(1)%patches(ipatch)
1016  bctype = patch%bcType
1017  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
1018  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
1019  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
1020  iregsrc = patch%srcRegion
1021  ipatchsrc = patch%srcPatch
1022  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
1023  gridsrc => regions(iregsrc)%levels(1)%grid
1024 
1025  IF (regions(iregsrc)%procid /= global%myProcid) THEN
1026  CALL rflo_exchangednoderecv( regions(ireg),regions(iregsrc), &
1027  patch,patchsrc,.false.,grid%xyz )
1028  CALL rflo_mgelliptedgesfra( regions(ireg),grid%boundMoved, &
1029  grid%allExternal,grid%edgeMoved, &
1030  grid%arcLen12,grid%arcLen34, &
1031  grid%arcLen56,gridold%xyzOld, &
1032  gridold%xyz,grid%xyz )
1033  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
1034  grid%edgeMoved,grid%arcLen12, &
1035  grid%arcLen34,grid%arcLen56, &
1036  gridold%xyzOld,grid%xyz )
1037  ENDIF
1038  ENDIF ! bcType
1039  ENDDO ! iPatch
1040 
1041  ENDIF ! region on this processor and active, grid moving
1042  ENDDO ! iReg
1043 
1044 ! - clear send requests
1045 
1046  DO ireg=1,global%nRegions
1047  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
1048  regions(ireg)%active==active .AND. & ! on my processor
1049  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
1050  CALL rflo_clearsendrequests( regions,ireg,.true. )
1051  ENDIF
1052  ENDDO
1053 
1054  ENDDO ! iPass
1055 
1056 ! finalize --------------------------------------------------------------------
1057 
1058  CALL deregisterfunction( global )
1059 
1060 END SUBROUTINE rflo_mgelliptinterfacesfra
1061 
1062 
1063 ! ******************************************************************************
1064 ! End
1065 ! ******************************************************************************
1066 
1067 END MODULE rflo_modmovegridelliptfra
1068 
1069 ! ******************************************************************************
1070 !
1071 ! RCS Revision history:
1072 !
1073 ! $Log: RFLO_ModMoveGridElliptFra.F90,v $
1074 ! Revision 1.16 2009/08/27 14:04:50 mtcampbe
1075 ! Updated to enable burning motion with symmetry boundaries and enhanced
1076 ! burnout code.
1077 !
1078 ! Revision 1.15 2008/12/06 08:44:16 mtcampbe
1079 ! Updated license.
1080 !
1081 ! Revision 1.14 2008/11/19 22:17:27 mtcampbe
1082 ! Added Illinois Open Source License/Copyright
1083 !
1084 ! Revision 1.13 2006/03/18 13:26:21 wasistho
1085 ! added orthDir and orthWghtX,Y,Z
1086 !
1087 ! Revision 1.12 2006/03/18 11:48:51 wasistho
1088 ! cosmetics
1089 !
1090 ! Revision 1.11 2006/03/18 11:04:10 wasistho
1091 ! screen printed global skewness and minvol
1092 !
1093 ! Revision 1.10 2006/03/18 09:24:47 wasistho
1094 ! modified mgElliptEdgesFra
1095 !
1096 ! Revision 1.9 2006/03/18 08:17:48 wasistho
1097 ! added mgElliptSurfacesFra
1098 !
1099 ! Revision 1.8 2006/03/16 08:27:09 wasistho
1100 ! added interfacesFra and edgesFra
1101 !
1102 ! Revision 1.7 2006/03/15 06:38:19 wasistho
1103 ! added region and global skewness
1104 !
1105 ! Revision 1.6 2006/03/14 04:41:13 wasistho
1106 ! added RFLO_EdgeDeformationStraight
1107 !
1108 ! Revision 1.5 2006/03/09 19:14:25 wasistho
1109 ! prepared for global elliptic grid motion
1110 !
1111 ! Revision 1.4 2006/03/05 19:08:31 wasistho
1112 ! set computational space coordinates from initial grid
1113 !
1114 ! Revision 1.3 2006/03/04 04:35:41 wasistho
1115 ! call elliptGridSmooRegion if blockShape is normal
1116 !
1117 ! Revision 1.2 2006/03/03 06:07:05 wasistho
1118 ! enabled global elliptic PDE grid motion
1119 !
1120 ! Revision 1.1 2006/02/11 03:42:14 wasistho
1121 ! added ModMoveGridElliptGlo/Fra
1122 !
1123 !
1124 !
1125 ! ******************************************************************************
1126 
1127 
1128 
1129 
1130 
1131 
1132 
1133 
1134 
1135 
1136 
subroutine rflo_edgedeformation(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine, public rflo_elliptgridsmoo(regions, resid)
subroutine rflo_copygeometrydummy(region)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jpnbeg
subroutine, public rflo_elliptgridsmooregion(region, resid)
subroutine rflo_calccellcentroids(region)
subroutine, public rflo_movegridelliptfra(regions)
double s
Definition: blastest.C:80
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)
subroutine rflo_mgelliptboundaries(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine rflo_exchangednoderecv(region, regionSrc, patch, patchSrc, average, dNode)
subroutine rflo_c2favgcoeffs(region)
double sqrt(double d)
Definition: double.h:73
subroutine rflo_changeinteriorgrid(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, xyz)
subroutine, public rflo_gridqualityglobal(regions)
subroutine, public rflo_mgelliptbnddeformation(region)
subroutine rflo_movegridinterfaces(regions)
subroutine rflo_calccontrolvolumes(region)
subroutine rflo_tfint1d(s, p1, p2, xyz)
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_edgedeformationstraight(region, boundMoved, edgeStraight, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOrig, xyzOld, dNode)
subroutine rflo_mgframeorthoshift(regions)
subroutine rflo_mgframemovecorners(regions)
Definition: patch.h:74
subroutine rflo_calcfacevectors(region)
subroutine rflo_exchangegeometry(regions)
subroutine rflo_generatecoarsegrids(region)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ipnbeg
void int int REAL * x
Definition: read.cpp:74
subroutine rflo_mgelliptsurfacesfra(regions, someMoved)
subroutine rflo_getdeformation(region, boundMoved, dNode)
subroutine rflo_clearsendrequests(regions, iReg, geometry)
subroutine rflo_calcgridspeeds(region)
subroutine rflo_mgelliptinterfacesfra(regions)
subroutine rflo_getdimensphysnodes(region, iLev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend)
subroutine rflo_exchangednodecopy(region, regionSrc, patch, patchSrc, average, dNode, dNodeSrc)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine rflo_mgelliptedgesfra(region, boundMoved, allExternal, edgeMoved, arcLen12, arcLen34, arcLen56, xyzRef, xyzOld, dNode)
subroutine rflo_boundarydeformation(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzOld, dNode)
subroutine grid(bp)
Definition: setup_py.f90:257
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)
RT a() const
Definition: Line_2.h:140
subroutine rflo_laplacegridsmoo(regions, resid)
subroutine rflo_calcfacecentroids(region)
subroutine rflo_mgframecorrectneighbors(regions)