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 | ! P5E.f90: set edge estimates via degree-5 polynomials. |
---|
31 | ! |
---|
32 | ! Darren Engwirda |
---|
33 | ! 25-Mar-2019 |
---|
34 | ! de2363 [at] columbia [dot] edu |
---|
35 | ! |
---|
36 | ! |
---|
37 | |
---|
38 | subroutine p5e(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,idof,head,tail |
---|
75 | logical :: okay |
---|
76 | real*8 :: xhat,fEPS |
---|
77 | real*8 :: delh(-3:+2) |
---|
78 | real*8 :: xmap(-3:+3) |
---|
79 | real*8 :: fhat(+6, nvar) |
---|
80 | real*8 :: ivec(+6,-3:+3) |
---|
81 | real*8 :: cmat(+6,+6) |
---|
82 | |
---|
83 | integer, parameter :: NSIZ = +6 |
---|
84 | real*8 , parameter :: ZERO = 1.d-14 |
---|
85 | |
---|
86 | head = +4 ; tail = npos - 3 |
---|
87 | |
---|
88 | if (npos.le.6) then |
---|
89 | !----- default to reduced order if insufficient points ! |
---|
90 | call p3e (npos,nvar,ndof, & |
---|
91 | & delx,fdat,bclo, & |
---|
92 | & bchi,edge,dfdx, & |
---|
93 | & opts,dmin) |
---|
94 | end if |
---|
95 | |
---|
96 | if (npos.le.6) return |
---|
97 | |
---|
98 | !------ impose value/slope B.C.'s about lower endpoint ! |
---|
99 | |
---|
100 | call pbc(npos,nvar,ndof,delx, & |
---|
101 | & fdat,bclo,edge,dfdx, & |
---|
102 | & -1 ,dmin) |
---|
103 | |
---|
104 | !------ impose value/slope B.C.'s about upper endpoint ! |
---|
105 | |
---|
106 | call pbc(npos,nvar,ndof,delx, & |
---|
107 | & fdat,bchi,edge,dfdx, & |
---|
108 | & +1 ,dmin) |
---|
109 | |
---|
110 | ! Reconstruct edge-centred 6th-order polynomials. Com- ! |
---|
111 | ! pute values/slopes at edges directly. Mid.-order ex- ! |
---|
112 | ! trapolation at endpoints. ! |
---|
113 | |
---|
114 | if (size(delx).eq.+1) then |
---|
115 | |
---|
116 | do ipos = head , tail |
---|
117 | |
---|
118 | !--------------- reconstruction: constant grid-spacing ! |
---|
119 | |
---|
120 | do ivar = 1, nvar |
---|
121 | |
---|
122 | edge(ivar,ipos) = & |
---|
123 | & + ( 1.d0 / 60.d0) * & |
---|
124 | & fdat(1,ivar,ipos-3) & |
---|
125 | & - ( 8.d0 / 60.d0) * & |
---|
126 | & fdat(1,ivar,ipos-2) & |
---|
127 | & + (37.d0 / 60.d0) * & |
---|
128 | & fdat(1,ivar,ipos-1) & |
---|
129 | & + (37.d0 / 60.d0) * & |
---|
130 | & fdat(1,ivar,ipos+0) & |
---|
131 | & - ( 8.d0 / 60.d0) * & |
---|
132 | & fdat(1,ivar,ipos+1) & |
---|
133 | & + ( 1.d0 / 60.d0) * & |
---|
134 | & fdat(1,ivar,ipos+2) |
---|
135 | |
---|
136 | dfdx(ivar,ipos) = & |
---|
137 | & - ( 1.d0 / 90.d0) * & |
---|
138 | & fdat(1,ivar,ipos-3) & |
---|
139 | & + ( 5.d0 / 36.d0) * & |
---|
140 | & fdat(1,ivar,ipos-2) & |
---|
141 | & - (49.d0 / 36.d0) * & |
---|
142 | & fdat(1,ivar,ipos-1) & |
---|
143 | & + (49.d0 / 36.d0) * & |
---|
144 | & fdat(1,ivar,ipos+0) & |
---|
145 | & - ( 5.d0 / 36.d0) * & |
---|
146 | & fdat(1,ivar,ipos+1) & |
---|
147 | & + ( 1.d0 / 90.d0) * & |
---|
148 | & fdat(1,ivar,ipos+2) |
---|
149 | |
---|
150 | dfdx(ivar,ipos) = & |
---|
151 | dfdx(ivar,ipos) / delx(+1) |
---|
152 | |
---|
153 | end do |
---|
154 | |
---|
155 | end do |
---|
156 | |
---|
157 | else |
---|
158 | |
---|
159 | fEPS = ZERO * dmin |
---|
160 | |
---|
161 | do ipos = head , tail |
---|
162 | |
---|
163 | !--------------- reconstruction: variable grid-spacing ! |
---|
164 | |
---|
165 | delh(-3) = & |
---|
166 | & max(delx(ipos-3),dmin) |
---|
167 | delh(-2) = & |
---|
168 | & max(delx(ipos-2),dmin) |
---|
169 | delh(-1) = & |
---|
170 | & max(delx(ipos-1),dmin) |
---|
171 | delh(+0) = & |
---|
172 | & max(delx(ipos+0),dmin) |
---|
173 | delh(+1) = & |
---|
174 | & max(delx(ipos+1),dmin) |
---|
175 | delh(+2) = & |
---|
176 | & max(delx(ipos+2),dmin) |
---|
177 | |
---|
178 | xhat = .5d0 * delh(-1) + & |
---|
179 | & .5d0 * delh(+0) |
---|
180 | |
---|
181 | xmap(-3) = -( delh(-3) & |
---|
182 | & + delh(-2) & |
---|
183 | & + delh(-1) ) / xhat |
---|
184 | xmap(-2) = -( delh(-2) & |
---|
185 | & + delh(-1) ) / xhat |
---|
186 | xmap(-1) = - delh(-1) / xhat |
---|
187 | xmap(+0) = + 0.d0 |
---|
188 | xmap(+1) = + delh(+0) / xhat |
---|
189 | xmap(+2) = +( delh(+0) & |
---|
190 | & + delh(+1) ) / xhat |
---|
191 | xmap(+3) = +( delh(+0) & |
---|
192 | & + delh(+1) & |
---|
193 | & + delh(+2) ) / xhat |
---|
194 | |
---|
195 | !--------------------------- calc. integral basis vec. ! |
---|
196 | |
---|
197 | do idof = -3, +3 |
---|
198 | |
---|
199 | ivec(1,idof) = & |
---|
200 | & xmap(idof) ** 1 / 1.0d+0 |
---|
201 | ivec(2,idof) = & |
---|
202 | & xmap(idof) ** 2 / 2.0d+0 |
---|
203 | ivec(3,idof) = & |
---|
204 | & xmap(idof) ** 3 / 3.0d+0 |
---|
205 | ivec(4,idof) = & |
---|
206 | & xmap(idof) ** 4 / 4.0d+0 |
---|
207 | ivec(5,idof) = & |
---|
208 | & xmap(idof) ** 5 / 5.0d+0 |
---|
209 | ivec(6,idof) = & |
---|
210 | & xmap(idof) ** 6 / 6.0d+0 |
---|
211 | |
---|
212 | end do |
---|
213 | |
---|
214 | !--------------------------- linear system: lhs matrix ! |
---|
215 | |
---|
216 | do idof = +1, +6 |
---|
217 | |
---|
218 | cmat(1,idof) = ivec(idof,-2) & |
---|
219 | & - ivec(idof,-3) |
---|
220 | cmat(2,idof) = ivec(idof,-1) & |
---|
221 | & - ivec(idof,-2) |
---|
222 | cmat(3,idof) = ivec(idof,+0) & |
---|
223 | & - ivec(idof,-1) |
---|
224 | cmat(4,idof) = ivec(idof,+1) & |
---|
225 | & - ivec(idof,+0) |
---|
226 | cmat(5,idof) = ivec(idof,+2) & |
---|
227 | & - ivec(idof,+1) |
---|
228 | cmat(6,idof) = ivec(idof,+3) & |
---|
229 | & - ivec(idof,+2) |
---|
230 | |
---|
231 | end do |
---|
232 | |
---|
233 | !--------------------------- linear system: rhs vector ! |
---|
234 | |
---|
235 | do ivar = +1, nvar |
---|
236 | |
---|
237 | fhat(+1,ivar) = & |
---|
238 | & delx(ipos-3) * & |
---|
239 | & fdat(+1,ivar,ipos-3) / xhat |
---|
240 | fhat(+2,ivar) = & |
---|
241 | & delx(ipos-2) * & |
---|
242 | & fdat(+1,ivar,ipos-2) / xhat |
---|
243 | fhat(+3,ivar) = & |
---|
244 | & delx(ipos-1) * & |
---|
245 | & fdat(+1,ivar,ipos-1) / xhat |
---|
246 | fhat(+4,ivar) = & |
---|
247 | & delx(ipos+0) * & |
---|
248 | & fdat(+1,ivar,ipos+0) / xhat |
---|
249 | fhat(+5,ivar) = & |
---|
250 | & delx(ipos+1) * & |
---|
251 | & fdat(+1,ivar,ipos+1) / xhat |
---|
252 | fhat(+6,ivar) = & |
---|
253 | & delx(ipos+2) * & |
---|
254 | & fdat(+1,ivar,ipos+2) / xhat |
---|
255 | |
---|
256 | end do |
---|
257 | |
---|
258 | !------------------------- factor/solve linear systems ! |
---|
259 | |
---|
260 | call slv_6x6(cmat,NSIZ,fhat, & |
---|
261 | & NSIZ,nvar,fEPS, & |
---|
262 | & okay) |
---|
263 | |
---|
264 | if (okay .eqv. .true.) then |
---|
265 | |
---|
266 | do ivar = +1, nvar |
---|
267 | |
---|
268 | edge(ivar,ipos) = fhat(1,ivar) |
---|
269 | |
---|
270 | dfdx(ivar,ipos) = fhat(2,ivar) & |
---|
271 | & / xhat |
---|
272 | |
---|
273 | end do |
---|
274 | |
---|
275 | else |
---|
276 | |
---|
277 | !------------------------- fallback if system singular ! |
---|
278 | |
---|
279 | # ifdef __PPR_WARNMAT__ |
---|
280 | |
---|
281 | write(*,*) & |
---|
282 | & "WARNING::P5E - matrix-is-singular!" |
---|
283 | |
---|
284 | # endif |
---|
285 | |
---|
286 | do ivar = +1, nvar |
---|
287 | |
---|
288 | edge(ivar,ipos) = & |
---|
289 | & fdat(1,ivar,ipos-1) * 0.5d+0 + & |
---|
290 | & fdat(1,ivar,ipos-0) * 0.5d+0 |
---|
291 | |
---|
292 | dfdx(ivar,ipos) = & |
---|
293 | & fdat(1,ivar,ipos-0) * 0.5d+0 - & |
---|
294 | & fdat(1,ivar,ipos-1) * 0.5d+0 |
---|
295 | |
---|
296 | dfdx(ivar,ipos) = & |
---|
297 | & dfdx(ivar,ipos) / xhat |
---|
298 | |
---|
299 | end do |
---|
300 | |
---|
301 | end if |
---|
302 | |
---|
303 | end do |
---|
304 | |
---|
305 | end if |
---|
306 | |
---|
307 | return |
---|
308 | |
---|
309 | end subroutine |
---|
310 | |
---|
311 | |
---|
312 | |
---|