Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLO_ModMoveGridElliptGlo.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 global grid-motion routines with elliptic PDE smoothings.
26 !
27 ! Description: The procedure is similar to MOVEGRID_GLOBAL method, except
28 ! the Laplacian smoothing is replaced by elliptic PDE smoothing.
29 !
30 ! Notes: None.
31 !
32 ! ******************************************************************************
33 !
34 ! $Id: RFLO_ModMoveGridElliptGlo.F90,v 1.18 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_movegridelliptglo, &
56 
57 ! private : RFLO_MgElliptSurfacesGlo
58 ! RFLO_MgElliptInterfacesGlo
59 ! RFLO_MgElliptEdgesGlo
60 
61 ! ******************************************************************************
62 ! Declarations and definitions
63 ! ******************************************************************************
64 
65  CHARACTER(CHRLEN) :: RCSIdentString = &
66  '$RCSfile: RFLO_ModMoveGridElliptGlo.F90,v $ $Revision: 1.18 $'
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 ! solving an elliptic PDE for node coordinates.
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_movegridelliptglo( regions )
93 
104 
105  IMPLICIT NONE
106 
107 #ifdef GENX
108  include 'roccomf90.h'
109 #endif
110 
111 ! ... parameters
112  TYPE(t_region), POINTER :: regions(:)
113 
114 ! ... loop variables
115  INTEGER :: ireg, iter, ipatch, ijk
116 
117 ! ... local variables
118  LOGICAL :: somemoved
119 
120  INTEGER :: bctype
121 
122  REAL(RFREAL) :: resid, globalresid
123  REAL(RFREAL), POINTER :: xyz(:,:), xyzold(:,:)
124 
125  TYPE(t_grid), POINTER :: grid, gridold
126  TYPE(t_global), POINTER :: global
127  TYPE(t_patch), POINTER :: patch
128 #ifdef GENX
129  DOUBLE PRECISION :: dalpha
130 #endif
131 
132 !******************************************************************************
133 
134  global => regions(1)%global
135 
136  CALL registerfunction( global,'RFLO_MoveGridElliptGlo',&
137  'RFLO_ModMoveGridElliptGlo.F90' )
138 
139 #ifdef GENX
140 ! update geometry buffers -----------------------------------------------------
141 
142  dalpha = global%dtMin/global%dTimeSystem
143  CALL com_call_function( global%genxHandleGm,1,dalpha )
144 #endif
145 
146 ! receive and distribute deformations for each region -------------------------
147 
148  CALL rflo_mgelliptsurfacesglo( regions,somemoved )
149 
150 ! fix interfaces between regions ----------------------------------------------
151 
152  IF (somemoved .eqv. .true.) THEN
153  CALL rflo_mgelliptinterfacesglo( regions )
154  ENDIF
155 
156 ! update grid, dummy, corner and edge cells -----------------------------------
157 
158  IF (somemoved) THEN
159  DO ireg=1,global%nRegions
160  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
161  regions(ireg)%active==active .AND. & ! on my processor
162  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
163 
164 ! ----- change the interior grid
165 
166  grid => regions(ireg)%levels(1)%grid
167  gridold => regions(ireg)%levels(1)%gridOld
168  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
169  grid%edgeMoved,grid%arcLen12, &
170  grid%arcLen34,grid%arcLen56, &
171  gridold%xyzOld,grid%xyz )
172 
173 ! ----- update coarse grids and dummy cells
174 
175  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
176  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
177  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
178  ENDIF ! region on this processor and active, grid moving
179  ENDDO ! iReg
180  CALL rflo_exchangegeometry( regions ) ! exchange geometry
181  ENDIF
182 
183 ! smooth grid by solving Laplace equation -------------------------------------
184 
185  IF (somemoved) THEN
186  resid = 0._rfreal
187  DO iter=1,global%moveGridNiter
188  CALL rflo_laplacegridsmoo( regions,resid )
189  ENDDO
190 ! DO iter=1,global%moveGridNiter
191 ! CALL RFLO_ElliptGridSmoo( regions,resid )
192 ! ENDDO
193 
194  IF (global%verbLevel /= verbose_none) THEN
195 #ifdef MPI
196  CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
197  masterproc,global%mpiComm,global%mpierr )
198  IF (global%mpierr /= 0) CALL errorstop( global,err_mpi_trouble,&
199  __line__ )
200 #else
201  globalresid = resid
202 #endif
203  IF (global%myProcid == masterproc .AND. &
204  global%verbLevel>=verbose_high) THEN
205  WRITE(stdout,2000) solver_name,global%skewness,global%minVol
206  WRITE(stdout,1000) solver_name,global%moveGridNiter,sqrt(globalresid)
207  ENDIF
208  ENDIF ! verbLevel
209  ENDIF ! someMoved
210 
211 ! update grid, dummy, corner and edge cells -----------------------------------
212 
213  DO ireg=1,global%nRegions
214  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
215  regions(ireg)%active==active .AND. & ! on my processor
216  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
217 
218 ! --- change xyz from coordinates to deformations
219 
220  xyz => regions(ireg)%levels(1)%grid%xyz
221  xyzold => regions(ireg)%levels(1)%gridOld%xyz
222 
223  DO ijk=lbound(xyz,2),ubound(xyz,2)
224  xyz(xcoord,ijk) = xyz(xcoord,ijk) - xyzold(xcoord,ijk)
225  xyz(ycoord,ijk) = xyz(ycoord,ijk) - xyzold(ycoord,ijk)
226  xyz(zcoord,ijk) = xyz(zcoord,ijk) - xyzold(zcoord,ijk)
227  ENDDO
228 
229 ! --- redistribute deformations at boundaries
230 
231  grid => regions(ireg)%levels(1)%grid
232  gridold => regions(ireg)%levels(1)%gridOld
233  grid%boundMoved(:) = .true.
234  grid%edgeMoved(:) = .true.
235 
236  DO ipatch=1,regions(ireg)%nPatches
237  patch => regions(ireg)%levels(1)%patches(ipatch)
238  bctype = patch%bcType
239 ! IF ((bcType>=BC_SYMMETRY .AND. bcType<=BC_SYMMETRY+BC_RANGE)) THEN
240 ! grid%boundMoved(patch%lbound) = .false.
241 ! ENDIF ! bcType
242  IF (bctype.eq.bc_symmetry) THEN
243  grid%boundMoved(patch%lbound) = .false.
244  ENDIF ! bcType
245  ENDDO ! iPatch
246 
247  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
248  grid%edgeMoved,grid%arcLen12, &
249  grid%arcLen34,grid%arcLen56, &
250  gridold%xyzOld,grid%xyz )
251 
252  CALL rflo_mgelliptbnddeformation( regions(ireg) )
253 
254 ! --- change xyz from deformations to coordinates
255 
256  CALL rflo_changeinteriorgrid( regions(ireg),grid%boundMoved, &
257  grid%edgeMoved,grid%arcLen12, &
258  grid%arcLen34,grid%arcLen56, &
259  gridold%xyzOld,grid%xyz )
260 
261 ! --- perform volume smoothing based on 3D Elliptic PDE
262 
263  IF (regions(ireg)%blockShape==region_shape_normal) THEN
264  DO iter=1,global%moveGridViter
265  CALL rflo_elliptgridsmooregion( regions(ireg),resid )
266  ENDDO
267  ENDIF
268 
269 ! --- update coarse grids and dummy cells
270 
271  CALL rflo_generatecoarsegrids( regions(ireg) ) ! coarsen finest grid
272  CALL rflo_copygeometrydummy( regions(ireg) ) ! copy to dummy nodes
273  CALL rflo_extrapolategeometry( regions(ireg) ) ! extrapolate
274  ENDIF ! region on this processor and active, grid moving
275  ENDDO ! iReg
276 
277  CALL rflo_exchangegeometry( regions ) ! exchange geometry
278 
279 ! print residual of 3D Elliptic PDE smoothing ---------------------------------
280 
281  IF (global%verbLevel /= verbose_none) THEN
282 #ifdef MPI
283  CALL mpi_reduce( resid,globalresid,1,mpi_rfreal,mpi_sum, &
284  masterproc,global%mpiComm,global%mpierr )
285  IF (global%mpierr /= 0) CALL errorstop( global,err_mpi_trouble,&
286  __line__ )
287 #else
288  globalresid = resid
289 #endif
290  IF (global%myProcid == masterproc .AND. &
291  global%verbLevel>=verbose_high) THEN
292  WRITE(stdout,1000) solver_name,global%moveGridViter,sqrt(globalresid)
293  ENDIF
294  ENDIF ! verbLevel
295 
296 ! calculate new metrics and grid speeds ---------------------------------------
297 
298  DO ireg=1,global%nRegions
299  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
300  regions(ireg)%active==active .AND. & ! on my processor
301  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
302  CALL rflo_calcfacevectors( regions(ireg) ) ! faces
303  CALL rflo_calccontrolvolumes( regions(ireg) ) ! volumes
304  CALL rflo_calccellcentroids( regions(ireg) ) ! cell 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  CALL rflo_calcgridspeeds( regions(ireg) ) ! grid speeds
310  ENDIF ! region on this processor and active, grid moving
311  ENDDO ! iReg
312 
313 ! global grid quality measure -------------------------------------------------
314 
315  CALL rflo_gridqualityglobal( regions )
316 
317 ! finalize --------------------------------------------------------------------
318 
319  CALL deregisterfunction( global )
320 
321 1000 FORMAT(a,1x,'Block-Elliptic-PDE grid motion:',i6,1pe13.4)
322 2000 FORMAT(a,1x,'global skewness, minvol:',2(1pe14.5))
323 
324 END SUBROUTINE rflo_movegridelliptglo
325 
326 
327 !******************************************************************************
328 !
329 ! Purpose: calculate node displacements on non-external flat patches
330 ! (finest grid only).
331 !
332 ! Description: none.
333 !
334 ! Input: region = grid dimensions
335 !
336 ! Output: xyz = updated deformations at boundaries.
337 !
338 !******************************************************************************
339 
340 SUBROUTINE rflo_mgelliptbnddeformation( region )
341 
347 
348  IMPLICIT NONE
349 
350 #include "Indexing.h"
351 
352 ! ... parameters
353 
354  TYPE(t_region) :: region
355 
356 ! ... loop variables
357  INTEGER :: ipatch, ijk
358 
359 ! ... local variables
360  REAL(RFREAL), POINTER :: xyz(:,:), xyzold(:,:)
361  TYPE(t_patch), POINTER :: patch
362  TYPE(t_global), POINTER :: global
363 
364 !******************************************************************************
365 
366  global => region%global
367 
368  CALL registerfunction( global,'RFLO_MgElliptBndDeformation',&
369  'RFLO_ModMoveGridElliptGlo.F90' )
370 
371 ! transform deformations to coordinates ---------------------------------------
372 
373  xyz => region%levels(1)%grid%xyz
374  xyzold => region%levels(1)%gridOld%xyz
375 
376  DO ijk=lbound(xyz,2),ubound(xyz,2)
377  xyz(xcoord,ijk) = xyz(xcoord,ijk) + xyzold(xcoord,ijk)
378  xyz(ycoord,ijk) = xyz(ycoord,ijk) + xyzold(ycoord,ijk)
379  xyz(zcoord,ijk) = xyz(zcoord,ijk) + xyzold(zcoord,ijk)
380  ENDDO
381 
382 ! move nodes on boundaries with active edges ----------------------------------
383 
384  DO ipatch=1,region%nPatches
385  patch => region%levels(1)%patches(ipatch)
386 
387  IF (patch%bcMotion/=bc_external) THEN
388 
389  IF (patch%bndFlat .EQV. .true.) THEN
390  CALL rflo_elliptgridjac2d( region,patch,ipatch )
391 ! CALL RFLO_ElliptGridGauss2D( region,patch,iPatch )
392 ! CALL RFLO_ElliptGridSOR2D( region,patch,iPatch )
393  ELSEIF ((.NOT. patch%bndFlat) .AND. (patch%dirFlat > 0)) THEN
394  CALL rflo_movegridqflatpatch( region,patch,ipatch )
395  ELSEIF ((.NOT. patch%bndFlat) .AND. (patch%dirFlat < 0)) THEN
396 !need test WRITE(STDOUT,*)'curved patch found, ireg, ipatch', &
397 !need test region%iRegionGlobal, iPatch
398 !need test WRITE(STDOUT,*)'try activate RFLO_MoveGridCurvedPatch'
399 !need test WRITE(STDOUT,*)'and turn off 2nd RFLO_BoundaryDeformation'
400 !need test CALL RFLO_MoveGridCurvedPatch( region,patch,iPatch )
401  ENDIF
402 
403  ENDIF ! not.external
404  ENDDO ! iPatch
405 
406 ! transform coordinates to deformations ---------------------------------------
407 
408  DO ijk=lbound(xyz,2),ubound(xyz,2)
409  xyz(xcoord,ijk) = xyz(xcoord,ijk) - xyzold(xcoord,ijk)
410  xyz(ycoord,ijk) = xyz(ycoord,ijk) - xyzold(ycoord,ijk)
411  xyz(zcoord,ijk) = xyz(zcoord,ijk) - xyzold(zcoord,ijk)
412  ENDDO
413 
414 ! finalize --------------------------------------------------------------------
415 
416  CALL deregisterfunction( global )
417 
418 END SUBROUTINE rflo_mgelliptbnddeformation
419 
420 
421 !******************************************************************************
422 !
423 ! Purpose: receive and distribute the deformations of surfaces
424 ! in block-wise manner.
425 !
426 ! Description: none.
427 !
428 ! Input: regions = data of all grid regions.
429 !
430 ! Output: regions%levels%grid%xyz = deformations at the boundaries
431 ! someMoved = parts of grid moved.
432 !
433 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
434 ! is applied to the finest grid first.
435 !
436 !******************************************************************************
437 
438 SUBROUTINE rflo_mgelliptsurfacesglo( regions,someMoved )
439 
440  USE moddatatypes
441  USE modglobal, ONLY : t_global
442  USE moddatastruct, ONLY : t_region
443  USE modgrid, ONLY : t_grid
447  USE moderror
448  USE modparameters
449  IMPLICIT NONE
450 
451 ! ... parameters
452  LOGICAL :: somemoved
453 
454  TYPE(t_region), POINTER :: regions(:)
455 
456 ! ... loop variables
457  INTEGER :: ireg
458 
459 ! ... local variables
460  TYPE(t_grid), POINTER :: grid, gridold
461  TYPE(t_global), POINTER :: global
462 
463 !******************************************************************************
464 
465  global => regions(1)%global
466 
467  CALL registerfunction( global,'RFLO_MgElliptSurfacesGlo',&
468  'RFLO_ModMoveGridElliptGlo.F90' )
469 
470 ! move grid separately for each region ----------------------------------------
471 
472  somemoved = .false.
473 
474  DO ireg=1,global%nRegions
475  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
476  regions(ireg)%active==active .AND. & ! on my processor
477  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
478 
479  grid => regions(ireg)%levels(1)%grid
480  gridold => regions(ireg)%levels(1)%gridOld
481  somemoved = .true.
482 
483 ! --- store the old grid
484 
485  gridold%indSvel = grid%indSvel
486  gridold%ipc = grid%ipc
487  gridold%jpc = grid%jpc
488  gridold%kpc = grid%kpc
489  gridold%xyz(:,:) = grid%xyz(:,:)
490  gridold%si(:,:) = grid%si(:,:)
491  gridold%sj(:,:) = grid%sj(:,:)
492  gridold%sk(:,:) = grid%sk(:,:)
493  gridold%vol(:) = grid%vol(:)
494 
495 ! --- calculate arclengths between boundaries
496 
497  CALL rflo_arclengthbounds( regions(ireg),gridold%xyzOld, &
498  grid%arcLen12,grid%arcLen34,grid%arcLen56 )
499 
500 ! --- get the boundary deformations
501 
502  CALL rflo_getdeformation( regions(ireg),grid%boundMoved,grid%xyz )
503 
504 ! --- calculate deformations at remaining edges
505 
506  grid%edgeMoved(:) = .false.
507  CALL rflo_mgelliptedgesglo( regions(ireg),grid%boundMoved, &
508  grid%edgeMoved, grid%arcLen12, &
509  grid%arcLen34, grid%arcLen56, &
510  gridold%xyzOld,gridold%xyz,grid%xyz )
511 
512 ! --- correct deformations at straight edges
513 
514  IF (global%moveGridNiter < 1) THEN
515  CALL rflo_edgedeformationstraight( regions(ireg),grid%boundMoved, &
516  grid%edgeStraight,grid%edgeMoved, &
517  grid%arcLen12,grid%arcLen34,grid%arcLen56, &
518  gridold%xyzOld,gridold%xyz,grid%xyz )
519  ENDIF
520 
521 ! --- calculate deformations at remaining boundaries
522 
523  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
524  grid%edgeMoved,grid%arcLen12, &
525  grid%arcLen34,grid%arcLen56, &
526  gridold%xyzOld,grid%xyz )
527 
528  ENDIF ! region on this processor and active, grid moving
529  ENDDO ! iReg
530 
531 ! finalize --------------------------------------------------------------------
532 
533  CALL deregisterfunction( global )
534 
535 END SUBROUTINE rflo_mgelliptsurfacesglo
536 
537 
538 !******************************************************************************
539 !
540 ! Purpose: exchange deformations between the regions as to ensure
541 ! matching grid nodes at the interfaces.
542 !
543 ! Description: none.
544 !
545 ! Input: regions = data of all grid regions, deformations.
546 !
547 ! Output: regions%levels%grid%xyz = deformations at the boundaries.
548 !
549 ! Notes: grid%xyz temporarily stores nodal displacements. The deformation
550 ! is applied to the finest grid first.
551 !
552 !******************************************************************************
553 
554 SUBROUTINE rflo_mgelliptinterfacesglo( regions )
555 
559  IMPLICIT NONE
560 
561 ! ... parameters
562  TYPE(t_region), POINTER :: regions(:)
563 
564 ! ... loop variables
565  INTEGER :: ireg, ipatch, ipass
566 
567 ! ... local variables
568  INTEGER :: bctype, iregsrc, ipatchsrc
569 
570  TYPE(t_grid), POINTER :: grid, gridold, gridsrc
571  TYPE(t_global), POINTER :: global
572  TYPE(t_patch), POINTER :: patch, patchsrc
573 
574 !******************************************************************************
575 
576  global => regions(1)%global
577 
578  CALL registerfunction( global,'RFLO_MgElliptInterfacesGlo',&
579  'RFLO_ModMoveGridElliptGlo.F90' )
580 
581 ! fix interfaces between regions ----------------------------------------------
582 
583  DO ipass=1,2
584 
585 ! - copy / send deformations
586 
587  DO ireg=1,global%nRegions
588  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
589  regions(ireg)%active==active .AND. & ! on my processor
590  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
591 
592  grid => regions(ireg)%levels(1)%grid
593  gridold => regions(ireg)%levels(1)%gridOld
594 
595  DO ipatch=1,regions(ireg)%nPatches
596  patch => regions(ireg)%levels(1)%patches(ipatch)
597  bctype = patch%bcType
598  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
599  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
600  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
601  iregsrc = patch%srcRegion
602  ipatchsrc = patch%srcPatch
603  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
604  gridsrc => regions(iregsrc)%levels(1)%grid
605 
606  IF (regions(iregsrc)%procid == global%myProcid) THEN
607  CALL rflo_exchangednodecopy( regions(ireg),regions(iregsrc), &
608  patch,patchsrc,.false., &
609  grid%xyz,gridsrc%xyz )
610  CALL rflo_mgelliptedgesglo( regions(ireg),grid%boundMoved, &
611  grid%edgeMoved, grid%arcLen12, &
612  grid%arcLen34, grid%arcLen56, &
613  gridold%xyzOld,gridold%xyz,grid%xyz )
614  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
615  grid%edgeMoved,grid%arcLen12, &
616  grid%arcLen34,grid%arcLen56, &
617  gridold%xyzOld,grid%xyz )
618  ELSE
619  CALL rflo_exchangednodesend( regions(ireg),regions(iregsrc), &
620  patch,grid%xyz )
621  ENDIF
622  ENDIF ! bcType
623  ENDDO ! iPatch
624 
625  ENDIF ! region on this processor and active, grid moving
626  ENDDO ! iReg
627 
628 ! - receive deformations
629 
630  DO ireg=1,global%nRegions
631  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
632  regions(ireg)%active==active .AND. & ! on my processor
633  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
634 
635  grid => regions(ireg)%levels(1)%grid
636  gridold => regions(ireg)%levels(1)%gridOld
637 
638  DO ipatch=1,regions(ireg)%nPatches
639  patch => regions(ireg)%levels(1)%patches(ipatch)
640  bctype = patch%bcType
641  IF ((bctype>=bc_regionconf .AND. bctype<=bc_regionconf+bc_range) .OR. &
642  (bctype>=bc_tra_peri .AND. bctype<=bc_tra_peri +bc_range) .OR. &
643  (bctype>=bc_rot_peri .AND. bctype<=bc_rot_peri +bc_range)) THEN
644  iregsrc = patch%srcRegion
645  ipatchsrc = patch%srcPatch
646  patchsrc => regions(iregsrc)%levels(1)%patches(ipatchsrc)
647  gridsrc => regions(iregsrc)%levels(1)%grid
648 
649  IF (regions(iregsrc)%procid /= global%myProcid) THEN
650  CALL rflo_exchangednoderecv( regions(ireg),regions(iregsrc), &
651  patch,patchsrc,.false.,grid%xyz )
652  CALL rflo_mgelliptedgesglo( regions(ireg),grid%boundMoved, &
653  grid%edgeMoved, grid%arcLen12, &
654  grid%arcLen34, grid%arcLen56, &
655  gridold%xyzOld,gridold%xyz,grid%xyz )
656  CALL rflo_boundarydeformation( regions(ireg),grid%boundMoved, &
657  grid%edgeMoved,grid%arcLen12, &
658  grid%arcLen34,grid%arcLen56, &
659  gridold%xyzOld,grid%xyz )
660  ENDIF
661  ENDIF ! bcType
662  ENDDO ! iPatch
663 
664  ENDIF ! region on this processor and active, grid moving
665  ENDDO ! iReg
666 
667 ! - clear send requests
668 
669  DO ireg=1,global%nRegions
670  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
671  regions(ireg)%active==active .AND. & ! on my processor
672  regions(ireg)%mixtInput%moveGrid) THEN ! and moving
673  CALL rflo_clearsendrequests( regions,ireg,.true. )
674  ENDIF
675  ENDDO
676 
677  ENDDO ! iPass
678 
679 ! finalize --------------------------------------------------------------------
680 
681  CALL deregisterfunction( global )
682 
683 END SUBROUTINE rflo_mgelliptinterfacesglo
684 
685 
686 !******************************************************************************
687 !
688 ! Purpose: calculate node displacements on those edges whose end points have
689 ! moved, but the associated boundaries were not updated yet (finest
690 ! grid only).
691 !
692 ! Description: points along an edge are shifted using 1-D linear transfinite
693 ! interpolation (TFI).
694 !
695 ! Input: region = grid dimensions
696 ! boundMoved = flag for boundaries of a region which have moved
697 ! arcLen12 = arclength between i=const. boundaries for each j, k
698 ! arcLen34 = arclength between j=const. boundaries for each k, i
699 ! arcLen56 = arclength between k=const. boundaries for each i, j
700 ! xyzOld = grid from previous time step.
701 !
702 ! Output: edgeMoved = flag if discretization at an edge was changed
703 ! dNode = updated deformations at edges.
704 !
705 ! Notes: variable dNode contains the whole 3-D field.
706 !
707 !******************************************************************************
708 
709 SUBROUTINE rflo_mgelliptedgesglo( region,boundMoved,edgeMoved, &
710  arclen12,arclen34,arclen56,xyzref,xyzold, &
711  dnode )
712 
713  USE moddatatypes
714  USE moddatastruct, ONLY : t_region
717  USE moderror
718  USE modparameters
719  IMPLICIT NONE
720 
721 #include "Indexing.h"
722 
723 ! ... parameters
724  LOGICAL :: boundmoved(6), edgemoved(12)
725 
726  REAL(RFREAL), POINTER :: arclen12(:,:), arclen34(:,:), arclen56(:,:)
727  REAL(RFREAL), POINTER :: dnode(:,:), xyzold(:,:), xyzref(:,:)
728 
729  TYPE(t_region) :: region
730 
731 ! ... loop variables
732  INTEGER :: iedge, ind
733 
734 ! ... local variables
735  INTEGER :: ilev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend, l1c, l2c
736  INTEGER :: indbeg, indend, ijkn, ijkn1, ijknbeg, ijknend, inoff, ijnoff
737  INTEGER :: switch(12,9)
738 
739  REAL(RFREAL) :: arclen, ds, s, dn(3), dnbeg(3), dnend(3)
740 
741 !******************************************************************************
742 
743  CALL registerfunction( region%global,'RFLO_MgElliptEdgesGlo',&
744  'RFLO_ModMoveGridElliptGlo.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 edge switch -------------------------------------------------------------
754 ! switch(:,1) = begins at boundary
755 ! switch(:,2) = ends on boundary
756 ! switch(:,3) = right boundary
757 ! switch(:,4) = left boundary
758 ! switch(:,5) = direction (from-to boundary)
759 ! switch(:,6) = start index
760 ! switch(:,7) = end index
761 ! switch(:,8) = constant index in 1st direction
762 ! switch(:,9) = constant index in 2nd direction
763 
764  switch( 1,:) = (/5, 6, 1, 3, 56, kpnbeg, kpnend, ipnbeg, jpnbeg/)
765  switch( 2,:) = (/3, 4, 1, 6, 34, jpnbeg, jpnend, kpnend, ipnbeg/)
766  switch( 3,:) = (/5, 6, 1, 4, 56, kpnbeg, kpnend, ipnbeg, jpnend/)
767  switch( 4,:) = (/3, 4, 1, 5, 34, jpnbeg, jpnend, kpnbeg, ipnbeg/)
768  switch( 5,:) = (/5, 6, 2, 3, 56, kpnbeg, kpnend, ipnend, jpnbeg/)
769  switch( 6,:) = (/3, 4, 2, 6, 34, jpnbeg, jpnend, kpnend, ipnend/)
770  switch( 7,:) = (/5, 6, 2, 4, 56, kpnbeg, kpnend, ipnend, jpnend/)
771  switch( 8,:) = (/3, 4, 2, 5, 34, jpnbeg, jpnend, kpnbeg, ipnend/)
772  switch( 9,:) = (/1, 2, 3, 5, 12, ipnbeg, ipnend, jpnbeg, kpnbeg/)
773  switch(10,:) = (/1, 2, 3, 6, 12, ipnbeg, ipnend, jpnbeg, kpnend/)
774  switch(11,:) = (/1, 2, 4, 5, 12, ipnbeg, ipnend, jpnend, kpnbeg/)
775  switch(12,:) = (/1, 2, 4, 6, 12, ipnbeg, ipnend, jpnend, kpnend/)
776 
777 ! edge movement flag ----------------------------------------------------------
778 
779  IF (boundmoved(1)) THEN
780  edgemoved( 1) = .true.; edgemoved( 2) = .true.
781  edgemoved( 3) = .true.; edgemoved( 4) = .true.
782  ENDIF
783  IF (boundmoved(2)) THEN
784  edgemoved( 5) = .true.; edgemoved( 6) = .true.
785  edgemoved( 7) = .true.; edgemoved( 8) = .true.
786  ENDIF
787  IF (boundmoved(3)) THEN
788  edgemoved( 1) = .true.; edgemoved( 5) = .true.
789  edgemoved( 9) = .true.; edgemoved(10) = .true.
790  ENDIF
791  IF (boundmoved(4)) THEN
792  edgemoved( 3) = .true.; edgemoved( 7) = .true.
793  edgemoved(11) = .true.; edgemoved(12) = .true.
794  ENDIF
795  IF (boundmoved(5)) THEN
796  edgemoved( 4) = .true.; edgemoved( 8) = .true.
797  edgemoved( 9) = .true.; edgemoved(11) = .true.
798  ENDIF
799  IF (boundmoved(6)) THEN
800  edgemoved( 2) = .true.; edgemoved( 6) = .true.
801  edgemoved(10) = .true.; edgemoved(12) = .true.
802  ENDIF
803 
804 ! loop over all 12 edges ------------------------------------------------------
805 
806  DO iedge=1,12
807  IF ((boundmoved(switch(iedge,1)) .OR. boundmoved(switch(iedge,2))) .AND. &
808  ((.NOT.boundmoved(switch(iedge,3))) .OR. &
809  (.NOT.boundmoved(switch(iedge,4)))) .AND. &
810  (.NOT.edgemoved(iedge))) THEN
811 
812  edgemoved(iedge) = .true.
813 
814  ds = 0._rfreal
815  indbeg = switch(iedge,6)
816  indend = switch(iedge,7)
817  l1c = switch(iedge,8)
818  l2c = switch(iedge,9)
819  DO ind=indbeg+1,indend-1
820  IF (switch(iedge,5) == 12) THEN
821  ijkn = indijk(ind ,l1c,l2c,inoff,ijnoff)
822  ijkn1 = indijk(ind-1 ,l1c,l2c,inoff,ijnoff)
823  ijknbeg = indijk(indbeg,l1c,l2c,inoff,ijnoff)
824  ijknend = indijk(indend,l1c,l2c,inoff,ijnoff)
825  arclen = arclen12(l1c,l2c)
826  ELSE IF (switch(iedge,5) == 34) THEN
827  ijkn = indijk(l2c,ind ,l1c,inoff,ijnoff)
828  ijkn1 = indijk(l2c,ind-1 ,l1c,inoff,ijnoff)
829  ijknbeg = indijk(l2c,indbeg,l1c,inoff,ijnoff)
830  ijknend = indijk(l2c,indend,l1c,inoff,ijnoff)
831  arclen = arclen34(l1c,l2c)
832  ELSE IF (switch(iedge,5) == 56) THEN
833  ijkn = indijk(l1c,l2c,ind ,inoff,ijnoff)
834  ijkn1 = indijk(l1c,l2c,ind-1 ,inoff,ijnoff)
835  ijknbeg = indijk(l1c,l2c,indbeg,inoff,ijnoff)
836  ijknend = indijk(l1c,l2c,indend,inoff,ijnoff)
837  arclen = arclen56(l1c,l2c)
838  ENDIF
839  dnbeg(:) = dnode(:,ijknbeg) + xyzold(:,ijknbeg) - xyzref(:,ijknbeg)
840  dnend(:) = dnode(:,ijknend) + xyzold(:,ijknend) - xyzref(:,ijknend)
841 
842  ds = ds + sqrt((xyzref(xcoord,ijkn)-xyzref(xcoord,ijkn1))**2 + &
843  (xyzref(ycoord,ijkn)-xyzref(ycoord,ijkn1))**2 + &
844  (xyzref(zcoord,ijkn)-xyzref(zcoord,ijkn1))**2)
845  s = ds/arclen
846 
847  CALL rflo_tfint1d( s,dnbeg,dnend,dn )
848  dnode(:,ijkn) = dn(:) + xyzref(:,ijkn) - xyzold(:,ijkn)
849  ENDDO ! i
850 
851  ENDIF ! boundMoved
852  ENDDO ! iEdge
853 
854 ! finalize --------------------------------------------------------------------
855 
856  CALL deregisterfunction( region%global )
857 
858 END SUBROUTINE rflo_mgelliptedgesglo
859 
860 
861 ! ******************************************************************************
862 ! End
863 ! ******************************************************************************
864 
865 END MODULE rflo_modmovegridelliptglo
866 
867 
868 ! ******************************************************************************
869 !
870 ! RCS Revision history:
871 !
872 ! $Log: RFLO_ModMoveGridElliptGlo.F90,v $
873 ! Revision 1.18 2009/08/27 14:04:50 mtcampbe
874 ! Updated to enable burning motion with symmetry boundaries and enhanced
875 ! burnout code.
876 !
877 ! Revision 1.17 2008/12/06 08:44:16 mtcampbe
878 ! Updated license.
879 !
880 ! Revision 1.16 2008/11/19 22:17:27 mtcampbe
881 ! Added Illinois Open Source License/Copyright
882 !
883 ! Revision 1.15 2006/05/19 03:14:49 wasistho
884 ! commented MovedGridCurvedPatch call
885 !
886 ! Revision 1.14 2006/03/18 11:49:28 wasistho
887 ! called GridQualityGlobal
888 !
889 ! Revision 1.13 2006/03/18 11:03:45 wasistho
890 ! screen printed global skewness and minvol
891 !
892 ! Revision 1.12 2006/03/18 09:25:01 wasistho
893 ! modified mgElliptEdgesGlo
894 !
895 ! Revision 1.11 2006/03/18 08:18:12 wasistho
896 ! added mgElliptSurfacesGlo
897 !
898 ! Revision 1.10 2006/03/17 06:38:17 wasistho
899 ! added call to non-flat-surface smoother
900 !
901 ! Revision 1.9 2006/03/16 08:26:49 wasistho
902 ! added interfacesGlo and edgesGlo
903 !
904 ! Revision 1.8 2006/03/08 06:41:01 wasistho
905 ! added moveGridViter
906 !
907 ! Revision 1.7 2006/03/05 19:08:24 wasistho
908 ! set computational space coordinates from initial grid
909 !
910 ! Revision 1.6 2006/03/04 04:35:36 wasistho
911 ! call elliptGridSmooRegion if blockShape is normal
912 !
913 ! Revision 1.5 2006/03/03 05:23:09 wasistho
914 ! made public MgElliptBndDeformation
915 !
916 ! Revision 1.4 2006/03/03 04:14:53 wasistho
917 ! activated elliptGridSmooRegion
918 !
919 ! Revision 1.3 2006/03/02 00:23:44 wasistho
920 ! prepared elliptic pde grid motion
921 !
922 ! Revision 1.2 2006/02/23 21:35:17 wasistho
923 ! initialzed resid
924 !
925 ! Revision 1.1 2006/02/11 03:42:14 wasistho
926 ! added ModMoveGridElliptGlo/Fra
927 !
928 !
929 !
930 ! ******************************************************************************
931 
932 
933 
934 
935 
936 
937 
938 
939 
940 
941 
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_elliptgridsor2d(region, patch, iPatch)
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
subroutine rflo_c2eavgcoeffs(region)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine, public rflo_elliptgridjac2d(region, patch, iPatch)
subroutine rflo_extrapolategeometry(region)
subroutine rflo_exchangednoderecv(region, regionSrc, patch, patchSrc, average, dNode)
subroutine rflo_c2favgcoeffs(region)
double sqrt(double d)
Definition: double.h:73
subroutine rflo_mgelliptedgesglo(region, boundMoved, edgeMoved, arcLen12, arcLen34, arcLen56, xyzRef, xyzOld, dNode)
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)
Definition: patch.h:74
subroutine rflo_calcfacevectors(region)
subroutine rflo_exchangegeometry(regions)
subroutine rflo_movegridsurfaces(regions, someMoved)
subroutine rflo_mgelliptinterfacesglo(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_getdeformation(region, boundMoved, dNode)
subroutine rflo_clearsendrequests(regions, iReg, geometry)
subroutine rflo_calcgridspeeds(region)
subroutine rflo_getdimensphysnodes(region, iLev, ipnbeg, ipnend, jpnbeg, jpnend, kpnbeg, kpnend)
subroutine rflo_exchangednodecopy(region, regionSrc, patch, patchSrc, average, dNode, dNodeSrc)
subroutine, public rflo_movegridqflatpatch(region, patch, iPatch)
subroutine, public rflo_movegridelliptglo(regions)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
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
subroutine, public rflo_movegridcurvedpatch(region, patch, iPatch)
**********************************************************************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, public rflo_elliptgridgauss2d(region, patch, iPatch)
subroutine rflo_checkmetrics(iReg, region)
subroutine rflo_mgelliptsurfacesglo(regions, someMoved)
RT a() const
Definition: Line_2.h:140
subroutine rflo_laplacegridsmoo(regions, resid)