1 | |
---|
2 | ! |
---|
3 | ! This program may be freely redistributed under the |
---|
4 | ! condition that the copyright notices (including this |
---|
5 | ! entire header) are not removed, and no compensation |
---|
6 | ! is received through use of the software. Private, |
---|
7 | ! research, and institutional use is free. You may |
---|
8 | ! distribute modified versions of this code UNDER THE |
---|
9 | ! CONDITION THAT THIS CODE AND ANY MODIFICATIONS MADE |
---|
10 | ! TO IT IN THE SAME FILE REMAIN UNDER COPYRIGHT OF THE |
---|
11 | ! ORIGINAL AUTHOR, BOTH SOURCE AND OBJECT CODE ARE |
---|
12 | ! MADE FREELY AVAILABLE WITHOUT CHARGE, AND CLEAR |
---|
13 | ! NOTICE IS GIVEN OF THE MODIFICATIONS. Distribution |
---|
14 | ! of this code as part of a commercial system is |
---|
15 | ! permissible ONLY BY DIRECT ARRANGEMENT WITH THE |
---|
16 | ! AUTHOR. (If you are not directly supplying this |
---|
17 | ! code to a customer, and you are instead telling them |
---|
18 | ! how they can obtain it for free, then you are not |
---|
19 | ! required to make any arrangement with me.) |
---|
20 | ! |
---|
21 | ! Disclaimer: Neither I nor: Columbia University, the |
---|
22 | ! National Aeronautics and Space Administration, nor |
---|
23 | ! the Massachusetts Institute of Technology warrant |
---|
24 | ! or certify this code in any way whatsoever. This |
---|
25 | ! code is provided "as-is" to be used at your own risk. |
---|
26 | ! |
---|
27 | ! |
---|
28 | |
---|
29 | ! |
---|
30 | ! P1E.f90: set edge estimates via degree-1 polynomials. |
---|
31 | ! |
---|
32 | ! Darren Engwirda |
---|
33 | ! 09-Sep-2016 |
---|
34 | ! de2363 [at] columbia [dot] edu |
---|
35 | ! |
---|
36 | ! |
---|
37 | |
---|
38 | subroutine p1e(npos,nvar,ndof,delx, & |
---|
39 | & fdat,bclo,bchi,edge, & |
---|
40 | & dfdx,opts,dmin) |
---|
41 | |
---|
42 | ! |
---|
43 | ! NPOS no. edges over grid. |
---|
44 | ! NVAR no. state variables. |
---|
45 | ! NDOF no. degrees-of-freedom per grid-cell. |
---|
46 | ! DELX grid-cell spacing array. LENGTH(DELX) == +1 if |
---|
47 | ! spacing is uniform . |
---|
48 | ! FDAT grid-cell moments array. FDAT is an array with |
---|
49 | ! SIZE = NDOF-by-NVAR-by-NPOS-1 . |
---|
50 | ! BCLO boundary condition at lower endpoint. |
---|
51 | ! BCHI boundary condition at upper endpoint. |
---|
52 | ! EDGE edge-centred interp. for function-value. EDGE |
---|
53 | ! is an array with SIZE = NVAR-by-NPOS . |
---|
54 | ! DFDX edge-centred interp. for 1st-derivative. DFDX |
---|
55 | ! is an array with SIZE = NVAR-by-NPOS . |
---|
56 | ! OPTS method parameters. See RCON-OPTS for details . |
---|
57 | ! DMIN min. grid-cell spacing thresh . |
---|
58 | ! |
---|
59 | |
---|
60 | implicit none |
---|
61 | |
---|
62 | !------------------------------------------- arguments ! |
---|
63 | integer, intent( in) :: npos,nvar,ndof |
---|
64 | real*8 , intent( in) :: delx(:) |
---|
65 | real*8 , intent( in) :: fdat(:,:,:) |
---|
66 | type (rcon_ends), intent(in) :: bclo(:) |
---|
67 | type (rcon_ends), intent(in) :: bchi(:) |
---|
68 | real*8 , intent(out) :: edge(:,:) |
---|
69 | real*8 , intent(out) :: dfdx(:,:) |
---|
70 | real*8 , intent( in) :: dmin |
---|
71 | class(rcon_opts), intent(in) :: opts |
---|
72 | |
---|
73 | !------------------------------------------- variables ! |
---|
74 | integer :: ipos,ivar,head,tail |
---|
75 | real*8 :: dd10 |
---|
76 | real*8 :: delh(-1:+0) |
---|
77 | |
---|
78 | head = +2; tail = npos-1 |
---|
79 | |
---|
80 | if (npos.lt.2) return |
---|
81 | if (npos.eq.2) then |
---|
82 | !----- default to reduced order if insufficient points ! |
---|
83 | do ivar = 1,nvar |
---|
84 | |
---|
85 | edge(ivar,1) = fdat(1,ivar,1) |
---|
86 | dfdx(ivar,1) = 0.d0 |
---|
87 | |
---|
88 | edge(ivar,2) = fdat(1,ivar,1) |
---|
89 | dfdx(ivar,2) = 0.d0 |
---|
90 | |
---|
91 | end do |
---|
92 | end if |
---|
93 | |
---|
94 | if (npos.le.2) return |
---|
95 | |
---|
96 | ! Reconstruct edge-centred 2nd-order polynomials. Com- ! |
---|
97 | ! pute values/slopes at edges directly. Full-order ex- ! |
---|
98 | ! trapolation at endpoints. |
---|
99 | |
---|
100 | if (size(delx).eq.+1) then |
---|
101 | |
---|
102 | do ipos = head , tail |
---|
103 | |
---|
104 | !--------------- reconstruction: constant grid-spacing ! |
---|
105 | |
---|
106 | dd10 = delx(+1) * 2.d0 |
---|
107 | |
---|
108 | do ivar = +1, nvar |
---|
109 | |
---|
110 | edge(ivar,ipos) = & |
---|
111 | & + delx(+1) * & |
---|
112 | & fdat(1,ivar,ipos-1) & |
---|
113 | & + delx(+1) * & |
---|
114 | & fdat(1,ivar,ipos+0) |
---|
115 | |
---|
116 | dfdx(ivar,ipos) = & |
---|
117 | & - 2.0d+0 * & |
---|
118 | & fdat(1,ivar,ipos-1) & |
---|
119 | & + 2.0d+0 * & |
---|
120 | & fdat(1,ivar,ipos+0) |
---|
121 | |
---|
122 | edge(ivar,ipos) = & |
---|
123 | & edge(ivar,ipos) / dd10 |
---|
124 | dfdx(ivar,ipos) = & |
---|
125 | & dfdx(ivar,ipos) / dd10 |
---|
126 | |
---|
127 | end do |
---|
128 | |
---|
129 | end do |
---|
130 | |
---|
131 | else |
---|
132 | |
---|
133 | do ipos = head , tail |
---|
134 | |
---|
135 | !--------------- reconstruction: variable grid-spacing ! |
---|
136 | |
---|
137 | delh(-1) = & |
---|
138 | & max(delx(ipos-1),dmin) |
---|
139 | delh(+0) = & |
---|
140 | & max(delx(ipos+0),dmin) |
---|
141 | |
---|
142 | dd10 = delh(-1)+delh(+0) |
---|
143 | |
---|
144 | do ivar = +1, nvar |
---|
145 | |
---|
146 | edge(ivar,ipos) = & |
---|
147 | & + delh(+0) * & |
---|
148 | & fdat(1,ivar,ipos-1) & |
---|
149 | & + delh(-1) * & |
---|
150 | & fdat(1,ivar,ipos+0) |
---|
151 | |
---|
152 | dfdx(ivar,ipos) = & |
---|
153 | & - 2.0d+0 * & |
---|
154 | & fdat(1,ivar,ipos-1) & |
---|
155 | & + 2.0d+0 * & |
---|
156 | & fdat(1,ivar,ipos+0) |
---|
157 | |
---|
158 | edge(ivar,ipos) = & |
---|
159 | & edge(ivar,ipos) / dd10 |
---|
160 | dfdx(ivar,ipos) = & |
---|
161 | & dfdx(ivar,ipos) / dd10 |
---|
162 | |
---|
163 | end do |
---|
164 | |
---|
165 | end do |
---|
166 | |
---|
167 | end if |
---|
168 | |
---|
169 | !------------- 1st-order value/slope BC's at endpoints ! |
---|
170 | |
---|
171 | do ivar = +1, nvar |
---|
172 | |
---|
173 | edge(ivar,head-1) = & |
---|
174 | & fdat(+1,ivar,head-1) |
---|
175 | edge(ivar,tail+1) = & |
---|
176 | & fdat(+1,ivar,tail+0) |
---|
177 | |
---|
178 | dfdx(ivar,head-1) = 0.d0 |
---|
179 | dfdx(ivar,tail+1) = 0.d0 |
---|
180 | |
---|
181 | end do |
---|
182 | |
---|
183 | return |
---|
184 | |
---|
185 | end subroutine |
---|
186 | |
---|
187 | |
---|
188 | |
---|