2D AMRCLAW
Functions/Subroutines
ginit.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine ginit (msave, first, nvar, naux, start_time)
 Initializes solution on all grids at level of grid msave, by calling qinit(). More...
 

Function/Subroutine Documentation

◆ ginit()

subroutine ginit (   msave,
logical  first,
  nvar,
  naux,
  start_time 
)

Initializes solution on all grids at level of grid msave, by calling qinit().

If first = true, (first call to init), then allocate the soln storage area too, else was already allocated.

Parameters
msaveFirst grid in the grid list that stores all grids at this level
firstIs this the first call to init?
nvarnumber of equations for the system
nauxnumber of auxiliary variables
start_timestart time of current simulation

Definition at line 18 of file ginit.f.

References amr_module::alloc, amr_module::cornxlo, amr_module::cornylo, amr_module::hxposs, amr_module::hyposs, igetsp(), amr_module::levelptr, amr_module::mxnest, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::needs_to_be_set, amr_module::nestlevel, amr_module::nghost, amr_module::node, qinit(), amr_module::rnode, setaux(), amr_module::store1, amr_module::store2, amr_module::storeaux, and amr_module::timemult.

Referenced by domain(), and setgrd().

18 c
19  use amr_module
20  implicit double precision (a-h,o-z)
21  logical first
22 
23 
24 c ::::::::::::::::::::::::::::: GINIT ::::::::::::::::::::::::
25 c
26 c initializes soln on all grids at 'level' by calling qinit
27 c if first = true, (first call to init), then allocate the
28 c soln storage area too, else was already allocated.
29 c
30 c :::::::::::::::::::::::::::::::::::::::;::::::::::::::::::::
31 
32  if (msave .eq. 0) go to 99
33 
34  level = node(nestlevel,msave)
35  hx = hxposs(level)
36  hy = hyposs(level)
37  mptr = msave
38 
39  10 nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
40  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
41  mitot = nx + 2*nghost
42  mjtot = ny + 2*nghost
43  corn1 = rnode(cornxlo,mptr)
44  corn2 = rnode(cornylo,mptr)
45  if(.not. (first)) go to 20
46  loc = igetsp(mitot*mjtot*nvar)
47  node(store1,mptr) = loc
48  if (naux .gt. 0) then
49  locaux = igetsp(mitot*mjtot*naux)
50  do k = 1, mitot*mjtot*naux,naux ! set first component of aux to signal that it
51  alloc(locaux+k-1) = needs_to_be_set ! needs val, wasnt copied from other grids
52  end do
53 
54  call setaux(nghost,nx,ny,corn1,corn2,hx,hy,
55  & naux,alloc(locaux))
56  else
57  locaux = 1
58  endif
59  node(storeaux,mptr) = locaux
60  if (level .lt. mxnest) then
61  loc2 = igetsp(mitot*mjtot*nvar)
62  node(store2,mptr) = loc2
63  endif
64  rnode(timemult, mptr) = start_time
65  go to 30
66  20 continue
67 c
68 c if 2nd time through, put initial values in store2 so finer grids
69 c can be advanced with interpolation of their boundary values.
70 c new time soln should still be in location store1.
71 c
72  loc = node(store2,mptr)
73  locaux = node(storeaux,mptr)
74 c
75  30 continue
76  call qinit(nvar,nghost,nx,ny,corn1,corn2,hx,hy,
77  & alloc(loc),naux,alloc(locaux))
78 
79 c
80  mptr = node(levelptr, mptr)
81  if (mptr .ne. 0) go to 10
82 c
83 c
84  99 continue
85  return
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
function igetsp(nwords)
Allocate contiguous space of length nword in main storage array alloc.
Definition: igetsp.f:9
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter nestlevel
AMR level of the grid.
Definition: amr_module.f90:44
subroutine qinit(meqn, mbc, mx, my, xlower, ylower, dx, dy, q, maux, aux)
Definition: qinit.f90:3
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer, parameter ndilo
global i index of left border of this grid
Definition: amr_module.f90:108
real(kind=8), parameter needs_to_be_set
Definition: amr_module.f90:168
integer, parameter ndjlo
global j index of lower border of this grid
Definition: amr_module.f90:114
integer, parameter store1
pointer to the address of memory storing the first copy of solution data on this grid, usually for storing new solution
Definition: amr_module.f90:101
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
Definition: amr_module.f90:105
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
integer, parameter cornylo
y-coordinate of the lower border of this grid
Definition: amr_module.f90:145
integer, parameter levelptr
node number (index) of next grid on the same level
Definition: amr_module.f90:35
integer mxnest
Definition: amr_module.f90:198
integer nghost
Definition: amr_module.f90:232
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)
Definition: setaux.f90:2
integer, parameter storeaux
pointer to the address of memory storing auxiliary data on this grid
Definition: amr_module.f90:120
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218
Here is the call graph for this function:
Here is the caller graph for this function: