Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RungeKuttaMP.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: calculate solution at a new time level.
26 !
27 ! Description: the governing equations are integrated in time using
28 ! the classical 4-stage Runge-Kutta method (4th-order in
29 ! time) in low-storage formulation.
30 !
31 ! Input: regions = data of all regions.
32 !
33 ! Output: regions%levels%mixt = new solution after one time step.
34 !
35 ! Notes: none.
36 !
37 ! ******************************************************************************
38 !
39 ! $Id: RungeKuttaMP.F90,v 1.20 2008/12/06 08:44:10 mtcampbe Exp $
40 !
41 ! Copyright: (c) 2003-2006 by the University of Illinois
42 !
43 ! ******************************************************************************
44 
45 SUBROUTINE rungekuttamp( regions )
46 
47  USE moddatatypes
48  USE moddatastruct, ONLY : t_region
49  USE modglobal, ONLY : t_global
50  USE moderror
51  USE modparameters
52 
53 #ifdef RFLU
59  USE rflu_modmpi
60  USE rflu_modnscbc
62  USE rflu_modtime, ONLY: rflu_settimerk
63 #endif
64 
65  USE modinterfaces, ONLY : afterupdatemp, &
72 #ifdef RFLU
77 #ifdef GENX
79 #endif
80 
81 #ifdef PLAG
83 #endif
84 #endif
85 
86  IMPLICIT NONE
87 
88 ! ... parameters
89  TYPE(t_region), POINTER :: regions(:)
90 
91 ! ... loop variables
92  INTEGER :: ireg, ireglocal, istage
93 
94 ! ... local variables
95  INTEGER :: flowmodel
96 
97  TYPE(t_global), POINTER :: global
98  TYPE(t_region), POINTER :: pregion
99 
100 !******************************************************************************
101 
102  global => regions(1)%global
103 
104  CALL registerfunction( global,'RungeKuttaMP',&
105  'RungeKuttaMP.F90' )
106 
107 ! loop over stages and regions ================================================
108 
109 
110  DO istage=1,global%nrkSteps
111 #ifdef RFLO
112  DO ireg=1,global%nRegions
113  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
114  regions(ireg)%active==active) THEN ! on my processor
115 #endif
116 #ifdef RFLU
117  DO ireglocal=1,global%nRegionsLocal
118  ireg = ireglocal
119 #endif
120 
121 ! ----- set pointer and get models --------------------------------------------
122 
123  pregion => regions(ireg)
124 
125  flowmodel = regions(ireg)%mixtInput%flowModel
126  regions(ireg)%irkStep = istage
127 
128 #ifdef RFLU
129 ! ----- Set RK time -----------------------------------------------------------
130 
131  CALL rflu_settimerk(pregion,istage)
132 
133 ! ----- RFLU fill GENX incoming buffers ---------------------------------------
134 
135 #ifdef GENX
136  CALL rflu_genx_initbflag(pregion)
137 #endif
138  CALL rflu_updateboundaryvalues(regions(ireg),istage)
139 
140 ! ----- Scale grid speeds -----------------------------------------------------
141 
142  CALL rflu_setgridspeedscalefactor(pregion)
143  CALL rflu_scalegridspeeds(pregion)
144 #endif
145 
146 ! ----- store previous solution; set dissipation to zero ----------------------
147 
148  CALL rkinitmp( regions(ireg),istage )
149 
150 ! ----- compute cell-gradients for higher-order scheme ------------------------
151 
152  CALL cellgradientsmp( regions(ireg) )
153 
154 ! ----- compute numerical dissipation -----------------------------------------
155 
156  CALL numericaldissipationmp( regions(ireg) )
157 
158 ! ----- compute viscous fluxes ------------------------------------------------
159 
160  IF ( flowmodel == flow_navst ) THEN
161  CALL viscousfluxesmp( regions(ireg) )
162  END IF ! flowModel
163 
164 ! ----- compute convective fluxes; form residual ------------------------------
165 
166  CALL convectivefluxesmp( regions(ireg) )
167 
168 ! ----- zero residuals --------------------------------------------------------
169 
170  CALL zeroresidualsmp(regions(ireg))
171 
172 ! ----- add source terms ------------------------------------------------------
173 
174  CALL sourcetermsmp( regions(ireg) )
175 
176 #ifdef RFLU
177 ! ----- add Equilibrium Eulerian corrections ----------------------------------
178 
179  CALL rflu_equilibriumeulerian( pregion )
180 #endif
181 
182 ! ----- zero out residuals in dummy cells -------------------------------------
183 
184  CALL zerodummycellsmp( regions(ireg) )
185 #ifdef RFLO
186  ENDIF ! region on this processor and active
187 #endif
188 #ifdef RFLU
189 ! ----- Descale grid speeds -----------------------------------------------------
190  CALL rflu_descalegridspeeds(pregion)
191 #endif
192  ENDDO ! iReg
193 
194 
195 #ifdef RFLU
196  IF(global%zoomFactor > 1) THEN
197  CALL rflu_timezoomdriver(regions)
198  ENDIF ! global%zoomFactor
199 #endif
200 
201 
202 #ifdef RFLO
203  DO ireg=1,global%nRegions
204  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
205  regions(ireg)%active==active) THEN ! on my processor
206 #endif
207 #ifdef RFLU
208  DO ireglocal=1,global%nRegionsLocal
209  ireg = ireglocal
210 #endif
211 
212 ! ----- set Region pointer
213 
214  pregion => regions(ireg)
215 
216 ! ----- update solution; sum up residuals -------------------------------------
217 
218  CALL rkupdatemp( regions(ireg),ireg,istage )
219 
220 #ifdef RFLU
221  IF ( rflu_nscbc_decidehavenscbc(pregion) .EQV. .true. ) THEN
222  CALL rflu_bxv_computevarscv(pregion)
223  END IF ! RFLU_NSCBC_DecideHaveNSCBC
224 #endif
225 
226 ! ----- perform checks and enforce after-update conditions --------------------
227 
228  CALL afterupdatemp( pregion,istage )
229 
230 
231 
232 #ifdef RFLO
233 ! ----- initiate communication kernel -----------------------------------------
234 
235  CALL initcommunicationmp( regions,ireg,istage )
236 #endif
237 
238 ! ----- update dependent variables --------------------------------------------
239 
240 #ifdef RFLU
241  CALL rflu_mpi_isendwrapper(pregion)
242  CALL rflu_setvarscontwrapper(pregion,1,pregion%grid%nCells)
243 
244 ! ----- update dependent variables on boundary faces --------------------------
245 
246  IF ( rflu_nscbc_decidehavenscbc(pregion) .EQV. .true. ) THEN
247  CALL rflu_bxv_setdependentvars(pregion)
248  END IF !
249 
250 #endif
251 #ifdef RFLO
252  CALL updatedependentvarsmp(regions(ireg))
253 #endif
254 
255 #ifdef RFLO
256  ENDIF ! region on this processor and active
257 #endif
258  ENDDO ! iReg
259 
260 #ifdef RFLO
261 ! - facilitate global non-dummy communication ---------------------------------
262 
263  CALL globalcommunicationmp( regions )
264 
265 ! - update boundary conditions ------------------------------------------------
266 
267  CALL updateboundaryconditionsmp( regions,istage )
268 #endif
269 
270 #ifdef RFLU
271  CALL rflu_mpi_copywrapper(regions)
272 
273  DO ireg = 1,global%nRegionsLocal
274  pregion => regions(ireg)
275 
276  CALL rflu_mpi_recvwrapper(pregion)
277  CALL rflu_setvarscontwrapper(pregion,pregion%grid%nCells+1, &
278  pregion%grid%nCellsTot)
279  CALL rflu_relp_transformwrapper(pregion)
280  END DO ! iReg
281 
282  DO ireg = 1,global%nRegionsLocal
283  pregion => regions(ireg)
284 
285  CALL rflu_mpi_clearrequestwrapper(pregion)
286  END DO ! iReg
287 
288 #ifdef PLAG
289  IF ( global%plagUsed .EQV. .true. ) THEN
290  CALL plag_rflu_commdriver(regions)
291 
292  DO ireg = 1,global%nRegionsLocal
293  pregion => regions(ireg)
294 
295  CALL rflu_setvarsdiscwrapper(pregion)
296  END DO ! iReg
297  END IF ! global%plagUsed
298 #endif
299 #endif
300  END DO ! istage
301 
302 ! finalize ====================================================================
303 
304  CALL deregisterfunction( global )
305 
306 END SUBROUTINE rungekuttamp
307 
308 ! ******************************************************************************
309 !
310 ! RCS Revision history:
311 !
312 ! $Log: RungeKuttaMP.F90,v $
313 ! Revision 1.20 2008/12/06 08:44:10 mtcampbe
314 ! Updated license.
315 !
316 ! Revision 1.19 2008/11/19 22:17:23 mtcampbe
317 ! Added Illinois Open Source License/Copyright
318 !
319 ! Revision 1.18 2007/04/20 16:07:48 mtcampbe
320 ! Updating for burnout support function RFLU_GENX_InitBFLAG
321 !
322 ! Revision 1.17 2007/04/14 14:23:24 mtcampbe
323 ! Updated for TZ
324 !
325 ! Revision 1.16 2007/03/27 00:39:50 haselbac
326 ! Removed call to PLAG_CalcnPclsTotGlobal, now in RFLU_TimeStepping
327 !
328 ! Revision 1.15 2007/03/20 22:02:29 fnajjar
329 ! Included call to PLAG_CalcnPclsTotGlobal
330 !
331 ! Revision 1.14 2006/08/21 16:10:01 haselbac
332 ! Adapted to name change
333 !
334 ! Revision 1.13 2006/08/19 15:48:25 mparmar
335 ! Added computations of boundary Cv and Dv for NSCBC implementation
336 !
337 ! Revision 1.12 2006/08/18 21:09:27 fnajjar
338 ! Removed IF around PLAG_RFLU_CommDriver for serial periodic cases
339 !
340 ! Revision 1.11 2006/03/25 21:40:03 haselbac
341 ! Added call to transforming data on related patches, cosmetics
342 !
343 ! Revision 1.10 2005/12/03 19:44:54 haselbac
344 ! Apparent bug fix: Separated call to RFLU_MPI_ClearRequestWrapper into separate loop
345 !
346 ! Revision 1.9 2005/12/01 21:52:14 fnajjar
347 ! Added IF statement around PLAG_RFLU_CommDriver, only active for more than one nRegions
348 !
349 ! Revision 1.8 2005/11/10 22:21:07 fnajjar
350 ! ACH: Proper fix for updating PLAG dv
351 !
352 ! Revision 1.7 2005/11/10 16:51:28 fnajjar
353 ! Added plagUsed IF statement around PLAG routines
354 !
355 ! Revision 1.6 2005/11/02 14:53:24 haselbac
356 ! Fady: Temporary fix so comm particles get non-cv vars updated properly
357 !
358 ! Revision 1.5 2005/05/18 22:04:41 fnajjar
359 ! Added PLAG communication routines; only initial implementation
360 !
361 ! Revision 1.4 2005/04/29 00:06:09 haselbac
362 ! Added routines to clear send requests
363 !
364 ! Revision 1.3 2005/04/15 15:06:06 haselbac
365 ! Converted to MPI
366 !
367 ! Revision 1.2 2005/03/31 16:31:02 haselbac
368 ! Added call to RFLU_SetTimeRK
369 !
370 ! Revision 1.1 2004/12/01 16:51:13 haselbac
371 ! Initial revision after changing case
372 !
373 ! Revision 1.18 2004/11/14 19:36:23 haselbac
374 ! Replaced call to UpdateDependentVarsMP by RFLU_SetVarsWrapper
375 !
376 ! Revision 1.17 2004/07/30 22:47:33 jferry
377 ! Implemented Equilibrium Eulerian method for Rocflu
378 !
379 ! Revision 1.16 2004/04/14 02:07:02 haselbac
380 ! Added grid-speed scaling calls for RFLU
381 !
382 ! Revision 1.15 2004/03/25 21:14:20 jferry
383 ! changed AfterUpdate to call most subroutines only after final RK stage
384 !
385 ! Revision 1.14 2004/03/02 21:47:28 jferry
386 ! Added After Update interactions
387 !
388 ! Revision 1.13 2004/02/26 21:11:58 wasistho
389 ! added globalCommunication
390 !
391 ! Revision 1.12 2004/02/26 21:01:46 haselbac
392 ! Enclosed updateBoundaryConditionsMP within ifdef RFLO
393 !
394 ! Revision 1.11 2004/01/29 22:52:47 haselbac
395 ! Added calls to RFLU_EnforceBoundsWrapper and updateDependentVarsMP
396 !
397 ! Revision 1.10 2003/12/04 03:23:06 haselbac
398 ! Added call to CellGradientsMP and validity check
399 !
400 ! Revision 1.9 2003/11/25 21:01:45 haselbac
401 ! Added calls to RFLU_UpdateDummyCells and ZeroResidualsMP
402 !
403 ! Revision 1.8 2003/11/20 16:40:35 mdbrandy
404 ! Backing out RocfluidMP changes from 11-17-03
405 !
406 ! Revision 1.5 2003/10/03 20:42:07 haselbac
407 ! Added Rocflu calls
408 !
409 ! Revision 1.4 2003/10/01 23:52:09 jblazek
410 ! Corrected bug in moving noslip wall BC and grid speeds.
411 !
412 ! Revision 1.3 2003/05/15 02:57:02 jblazek
413 ! Inlined index function.
414 !
415 ! Revision 1.2 2003/04/10 01:22:41 jblazek
416 ! Got rid of pRegion in ViscousFluxesMP.
417 !
418 ! Revision 1.1 2003/03/28 19:42:55 fnajjar
419 ! Initial import for RocfluidMP
420 !
421 ! ******************************************************************************
422 
423 
424 
425 
426 
427 
428 
subroutine convectivefluxesmp(region)
subroutine zeroresidualsmp(region)
subroutine viscousfluxesmp(region)
subroutine, public rflu_setgridspeedscalefactor(pRegion)
subroutine initcommunicationmp(regions, iReg, istage)
subroutine, public rflu_mpi_isendwrapper(pRegion)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine cellgradientsmp(region)
subroutine, public rflu_mpi_clearrequestwrapper(pRegion)
subroutine sourcetermsmp(region)
subroutine, public rflu_settimerk(pRegion, iStage)
LOGICAL function, public rflu_nscbc_decidehavenscbc(pRegion)
subroutine updateboundaryconditionsmp(regions, istage)
subroutine, public rflu_genx_initbflag(pRegion)
subroutine rungekuttamp(regions)
subroutine, public rflu_mpi_recvwrapper(pRegion)
subroutine, public rflu_relp_transformwrapper(pRegion)
subroutine, public plag_rflu_commdriver(regions)
subroutine zerodummycellsmp(region)
subroutine, public rflu_bxv_setdependentvars(pRegion)
subroutine rflu_setvarscontwrapper(pRegion, icgBeg, icgEnd)
subroutine rflu_updateboundaryvalues(region, istage)
subroutine rflu_equilibriumeulerian(pRegion)
subroutine globalcommunicationmp(regions)
subroutine rkupdatemp(region, iReg, istage)
Definition: RkUpdateMP.F90:45
subroutine, public rflu_descalegridspeeds(pRegion)
subroutine rflu_setvarsdiscwrapper(pRegion)
subroutine, public rflu_mpi_copywrapper(regions)
subroutine afterupdatemp(pRegion, istage)
subroutine, public rflu_bxv_computevarscv(pRegion)
subroutine, public rflu_scalegridspeeds(pRegion)
subroutine numericaldissipationmp(region)
subroutine rkinitmp(region, istage)
Definition: RkInitMP.F90:44
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine updatedependentvarsmp(region)
subroutine, public rflu_timezoomdriver(regions)