2D AMRCLAW
stst1.f
Go to the documentation of this file.
1 c
2 c --------------------------------------------------------------
3 c
4  subroutine stst1
5 c
6  use amr_module
7  implicit double precision (a-h,o-z)
8 
9 
10 c
11 c :::::::::::::::::::::::::::::: STST1 ::::::::::::::::::::::::::::::::
12 c intialize a few variables needed before calling user set up
13 c routine domain.
14 c the spatial and temporal stepsizes are set. the node array
15 c is kept as a linked list of free nodes. "ndfree" points to the
16 c head of the list, i.e. first free node. use first row of each
17 c col to hold this pointer, set by the macro "nextfree".
18 c the free space list, managed in lfree, will have first and
19 c last positions filled with an allocation of zero words,
20 c to avoid boundary cases.
21 c ::::::::::::::::::::::::::::::::::::::;::::::::::::::::::::::::::::::
22 c
23  ndfree = 1
24  do 10 i = 1, maxgr
25  node(nextfree,i) = i+1
26  10 continue
27 c
28 c the last free node will have a null pointer
29 
31 c
32 
33 c Initialize dynamic memory
34  call init_alloc()
35 
36  lfine = 1
37 C do 20 i = 1, memsize
38 C alloc(i) = WEIRD
39 C 20 continue
40 c
41 c initialize linked list of alloc storage as well.
42 c first and last locations are dummy placeholders of zero words
43 c of allocation each, to avoid boundary cases.
44 c
45  do 40 i = 1, lfdim
46  lfree(i,1) = 0
47  lfree(i,2) = 0
48  40 continue
49  lfree(3,1) =memsize + 2
50  lfree(2,1) = 1
51  lfree(2,2) =memsize
52  lenf = 3
53 
54 
55 c need to manage the boundary List too
56 c do i = 1, bndListSize
57 c bndList(i,nextfree) = i+1
58 c end do
59 c bndList(bndListSize,nextfree) = null
60 c ndfree_bnd = 1
61  call initbndrylist()
62 c
63 c after kcheck integrations of parent grid, move its refinements.
64 c finest level grid never needs to have its finer subgrids moved.
65 c
66  call inittimers() ! used to be done here, but needs to be called from restarting too when stst1 not called
67 
68  do 60 i = 1, maxlv
69  tvoll(i) = 0.d0
70  iregridcount(i) = 0
71  avenumgrids(i) = 0
72  numgrids(i) = 0
73  numcells(i) = 0
74  lstart(i) = 0
75  60 icheck(i) = 0
76 c
77 c finish initializing spatial and counting arrays
78 c
79  level = 2
80  70 if (level .gt. mxnest) go to 80
81  hxposs(level) = hxposs(level-1) / dble(intratx(level-1))
82  hyposs(level) = hyposs(level-1) / dble(intraty(level-1))
83  level = level + 1
84  go to 70
85  80 continue
86 
87 
88  return
89  end
90 c
91 c -------------------------------------------------------------------------
92 c
93  subroutine inittimers()
94 
95  use amr_module
96  !implicit double precision (a-h,o-z)
97 
98  timeflagger = 0
99  timebufnst = 0
100  timestepgrid = 0
101  timestepgridcpu = 0.d0
102  timebound = 0
103  timeboundcpu = 0.d0
104  timegrdfitall = 0
105  timeflagger = 0
106  timeflglvl = 0
107  timeflglvltot = 0
108  timegrdfit2 = 0
109  timeregridding = 0
110  timeregriddingcpu = 0.d0
111  timetick = 0
112  timetickcpu = 0.d0
113  timeupdating = 0
114  timevalout = 0
115  timevaloutcpu = 0.d0
116 
117  return
118  end
integer timetick
Definition: amr_module.f90:242
integer timegrdfit2
Definition: amr_module.f90:240
integer, dimension(maxlv) numcells
Definition: amr_module.f90:198
integer, parameter maxlv
Definition: amr_module.f90:174
integer ndfree
Definition: amr_module.f90:198
integer timegrdfitall
Definition: amr_module.f90:240
subroutine stst1
Definition: stst1.f:5
integer, parameter maxgr
Definition: amr_module.f90:173
integer timebound
Definition: amr_module.f90:241
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer timeflagger
Definition: amr_module.f90:242
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter lfdim
Definition: amr_module.f90:224
subroutine init_alloc()
Definition: init_alloc.f90:12
integer memsize
Definition: amr_module.f90:219
integer, dimension(maxlv) numgrids
Definition: amr_module.f90:198
integer, dimension(maxlv) iregridcount
Definition: amr_module.f90:238
subroutine inittimers()
Definition: stst1.f:94
real(kind=8) timevaloutcpu
Definition: amr_module.f90:245
integer timeupdating
Definition: amr_module.f90:239
real(kind=8), dimension(maxlv) avenumgrids
Definition: amr_module.f90:237
integer, dimension(maxlv) icheck
Definition: amr_module.f90:198
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
real(kind=8) timetickcpu
Definition: amr_module.f90:243
real(kind=8) timeboundcpu
Definition: amr_module.f90:244
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
integer timebufnst
Definition: amr_module.f90:242
real(kind=8) timeregriddingcpu
Definition: amr_module.f90:244
integer, dimension(lfdim, 2) lfree
Definition: amr_module.f90:225
integer, parameter nextfree
Definition: amr_module.f90:154
integer timestepgrid
Definition: amr_module.f90:241
integer lenf
Definition: amr_module.f90:225
real(kind=8) timestepgridcpu
Definition: amr_module.f90:244
integer, dimension(maxlv) intratx
Definition: amr_module.f90:198
integer, parameter null
Definition: amr_module.f90:155
integer mxnest
Definition: amr_module.f90:198
integer timeflglvl
Definition: amr_module.f90:240
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer timevalout
Definition: amr_module.f90:239
subroutine initbndrylist()
Definition: nodget.f:148
integer lfine
Definition: amr_module.f90:198
integer, dimension(maxlv) tvoll
Definition: amr_module.f90:238
integer timeregridding
Definition: amr_module.f90:239