1 | ! |
---|
2 | ! $Id$ |
---|
3 | ! |
---|
4 | C AGRIF (Adaptive Grid Refinement In Fortran) |
---|
5 | C |
---|
6 | C Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
7 | C Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
8 | C |
---|
9 | C This program is free software; you can redistribute it and/or modify |
---|
10 | C it under the terms of the GNU General Public License as published by |
---|
11 | C the Free Software Foundation; either version 2 of the License, or |
---|
12 | C (at your option) any later version. |
---|
13 | C |
---|
14 | C This program is distributed in the hope that it will be useful, |
---|
15 | C but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
16 | C MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
17 | C GNU General Public License for more details. |
---|
18 | C |
---|
19 | C You should have received a copy of the GNU General Public License |
---|
20 | C along with this program; if not, write to the Free Software |
---|
21 | C Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA. |
---|
22 | C |
---|
23 | C |
---|
24 | C |
---|
25 | CCC Module Agrif_Interpbasic |
---|
26 | C |
---|
27 | Module Agrif_Interpbasic |
---|
28 | C |
---|
29 | CCC Description: |
---|
30 | CCC Module containing different procedures of interpolation (linear,lagrange, |
---|
31 | CCC spline,...) used in the Agrif_Interpolation module. |
---|
32 | C |
---|
33 | C Modules used: |
---|
34 | USE Agrif_types |
---|
35 | C |
---|
36 | IMPLICIT NONE |
---|
37 | C |
---|
38 | Real,Dimension(Agrif_MaxRaff) :: tabdiff2, tabdiff3 |
---|
39 | Real,Dimension(5,Agrif_MaxRaff,3) :: tabppm |
---|
40 | Real,Dimension(:),Allocatable::tabtest4 |
---|
41 | Integer, Dimension(:,:), Allocatable :: indparent |
---|
42 | Integer, Dimension(:,:), Allocatable :: |
---|
43 | & indparentppm,indchildppm |
---|
44 | Integer, Dimension(:), Allocatable :: |
---|
45 | & indparentppm_1d,indchildppm_1d |
---|
46 | Real, Dimension(:,:),Allocatable :: coeffparent |
---|
47 | |
---|
48 | CONTAINS |
---|
49 | C Define procedures contained in this module |
---|
50 | C |
---|
51 | C ************************************************************************** |
---|
52 | CCC Subroutine Linear1d |
---|
53 | C ************************************************************************** |
---|
54 | C |
---|
55 | Subroutine Linear1d(x,y,np,nc, |
---|
56 | & s_parent,s_child,ds_parent,ds_child) |
---|
57 | C |
---|
58 | CCC Description: |
---|
59 | CCC Subroutine to do a linear 1D interpolation on a child grid (vector y) from |
---|
60 | CCC its parent grid (vector x). |
---|
61 | C |
---|
62 | CC Method: |
---|
63 | C |
---|
64 | C Declarations: |
---|
65 | C |
---|
66 | |
---|
67 | C |
---|
68 | C Arguments |
---|
69 | INTEGER :: np,nc |
---|
70 | REAL,INTENT(IN), DIMENSION(np) :: x |
---|
71 | REAL,INTENT(OUT), DIMENSION(nc) :: y |
---|
72 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
73 | C |
---|
74 | C Local scalars |
---|
75 | INTEGER :: i,coeffraf,locind_parent_left |
---|
76 | REAL :: ypos,globind_parent_left,globind_parent_right |
---|
77 | REAL :: invds, invds2 |
---|
78 | REAL :: ypos2,diff |
---|
79 | C |
---|
80 | C |
---|
81 | |
---|
82 | coeffraf = nint(ds_parent/ds_child) |
---|
83 | C |
---|
84 | if (coeffraf == 1) then |
---|
85 | C |
---|
86 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
87 | C |
---|
88 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
89 | C |
---|
90 | return |
---|
91 | C |
---|
92 | endif |
---|
93 | C |
---|
94 | ypos = s_child |
---|
95 | |
---|
96 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
97 | |
---|
98 | globind_parent_left = s_parent |
---|
99 | & + (locind_parent_left - 1)*ds_parent |
---|
100 | |
---|
101 | globind_parent_right = globind_parent_left + ds_parent |
---|
102 | |
---|
103 | C |
---|
104 | invds = 1./ds_parent |
---|
105 | invds2 = ds_child/ds_parent |
---|
106 | |
---|
107 | ypos2 = ypos*invds |
---|
108 | globind_parent_right=globind_parent_right*invds |
---|
109 | |
---|
110 | do i = 1,nc-1 |
---|
111 | C |
---|
112 | if (ypos2 > globind_parent_right) then |
---|
113 | locind_parent_left = locind_parent_left + 1. |
---|
114 | globind_parent_right = globind_parent_right + 1. |
---|
115 | endif |
---|
116 | |
---|
117 | diff=(globind_parent_right - ypos2) |
---|
118 | |
---|
119 | y(i) = (diff*x(locind_parent_left) |
---|
120 | & + (1.-diff)*x(locind_parent_left+1)) |
---|
121 | C |
---|
122 | ypos2 = ypos2 + invds2 |
---|
123 | C |
---|
124 | enddo |
---|
125 | C |
---|
126 | ypos = s_child + (nc-1)*ds_child |
---|
127 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
128 | C |
---|
129 | if (locind_parent_left == np) then |
---|
130 | C |
---|
131 | y(nc) = x(np) |
---|
132 | C |
---|
133 | else |
---|
134 | C |
---|
135 | globind_parent_left = s_parent |
---|
136 | & + (locind_parent_left - 1)*ds_parent |
---|
137 | C |
---|
138 | y(nc) = ((globind_parent_left + ds_parent - ypos) |
---|
139 | & *x(locind_parent_left) |
---|
140 | & + (ypos - globind_parent_left) |
---|
141 | & *x(locind_parent_left+1))*invds |
---|
142 | C |
---|
143 | endif |
---|
144 | C |
---|
145 | Return |
---|
146 | C |
---|
147 | C |
---|
148 | End Subroutine Linear1d |
---|
149 | |
---|
150 | Subroutine Linear1dprecompute2d(np2, np,nc, |
---|
151 | & s_parent,s_child,ds_parent,ds_child,dir) |
---|
152 | C |
---|
153 | CCC Description: |
---|
154 | CCC Subroutine to compute 2D coefficient and index for a linear 1D interpolation |
---|
155 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
156 | C |
---|
157 | CC Method: |
---|
158 | C |
---|
159 | C Declarations: |
---|
160 | C |
---|
161 | |
---|
162 | C |
---|
163 | C Arguments |
---|
164 | INTEGER :: np,nc,np2,dir |
---|
165 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
166 | C |
---|
167 | C Local scalars |
---|
168 | INTEGER :: i,j,coeffraf,locind_parent_left,inc,inc1,inc2,toto |
---|
169 | Integer, Dimension(:,:), Allocatable :: indparent_tmp |
---|
170 | Real, Dimension(:,:), Allocatable :: coeffparent_tmp |
---|
171 | REAL :: ypos,globind_parent_left,globind_parent_right |
---|
172 | REAL :: invds, invds2 |
---|
173 | REAL :: ypos2,diff |
---|
174 | C |
---|
175 | C |
---|
176 | |
---|
177 | coeffraf = nint(ds_parent/ds_child) |
---|
178 | C |
---|
179 | ypos = s_child |
---|
180 | |
---|
181 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
182 | |
---|
183 | globind_parent_left = s_parent |
---|
184 | & + (locind_parent_left - 1)*ds_parent |
---|
185 | |
---|
186 | globind_parent_right = globind_parent_left + ds_parent |
---|
187 | |
---|
188 | C |
---|
189 | invds = 1./ds_parent |
---|
190 | invds2 = ds_child/ds_parent |
---|
191 | |
---|
192 | ypos2 = ypos*invds |
---|
193 | globind_parent_right=globind_parent_right*invds |
---|
194 | |
---|
195 | if (.not.allocated(indparent)) then |
---|
196 | allocate(indparent(nc*np2,3),coeffparent(nc*np2,3)) |
---|
197 | else |
---|
198 | if (size(indparent,1)<nc*np2) then |
---|
199 | allocate(coeffparent_tmp(size(indparent,1),size(indparent,2))) |
---|
200 | allocate(indparent_tmp(size(indparent,1),size(indparent,2))) |
---|
201 | coeffparent_tmp=coeffparent |
---|
202 | indparent_tmp=indparent |
---|
203 | deallocate(indparent,coeffparent) |
---|
204 | allocate(indparent(nc*np2,3),coeffparent(nc*np2,3)) |
---|
205 | coeffparent(1:size(coeffparent_tmp,1),1:size(coeffparent_tmp,2)) |
---|
206 | & = coeffparent_tmp |
---|
207 | indparent(1:size(indparent_tmp,1),1:size(indparent_tmp,2)) |
---|
208 | & = indparent_tmp |
---|
209 | deallocate(indparent_tmp,coeffparent_tmp) |
---|
210 | endif |
---|
211 | endif |
---|
212 | |
---|
213 | do i = 1,nc-1 |
---|
214 | C |
---|
215 | if (ypos2 > globind_parent_right) then |
---|
216 | locind_parent_left = locind_parent_left + 1 |
---|
217 | globind_parent_right = globind_parent_right + 1. |
---|
218 | endif |
---|
219 | |
---|
220 | diff=(globind_parent_right - ypos2) |
---|
221 | indparent(i,dir) = locind_parent_left |
---|
222 | coeffparent(i,dir) = diff |
---|
223 | |
---|
224 | ypos2 = ypos2 + invds2 |
---|
225 | C |
---|
226 | enddo |
---|
227 | C |
---|
228 | ypos = s_child + (nc-1)*ds_child |
---|
229 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
230 | |
---|
231 | if (locind_parent_left == np) then |
---|
232 | indparent(nc,dir) = locind_parent_left-1 |
---|
233 | coeffparent(nc,dir) = 0. |
---|
234 | else |
---|
235 | C |
---|
236 | globind_parent_left = s_parent |
---|
237 | & + (locind_parent_left - 1)*ds_parent |
---|
238 | C |
---|
239 | indparent(nc,dir) = locind_parent_left |
---|
240 | |
---|
241 | coeffparent(nc,dir) = (globind_parent_left + ds_parent - ypos) |
---|
242 | & * invds |
---|
243 | endif |
---|
244 | |
---|
245 | do i=2, np2 |
---|
246 | inc = i*nc |
---|
247 | inc1 = (i-1)*nc |
---|
248 | inc2 = (i-2)*nc |
---|
249 | !CDIR ALTCODE |
---|
250 | indparent(1+inc1:inc,dir) = indparent(1+inc2:inc1,dir)+np |
---|
251 | !CDIR ALTCODE |
---|
252 | coeffparent(1+inc1:inc,dir) =coeffparent(1:nc,dir) |
---|
253 | enddo |
---|
254 | |
---|
255 | Return |
---|
256 | C |
---|
257 | C |
---|
258 | End Subroutine Linear1dprecompute2d |
---|
259 | |
---|
260 | |
---|
261 | |
---|
262 | Subroutine Linear1dprecompute(np,nc, |
---|
263 | & s_parent,s_child,ds_parent,ds_child) |
---|
264 | C |
---|
265 | CCC Description: |
---|
266 | CCC Subroutine to compute 1D coefficient and index for a linear 1D interpolation |
---|
267 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
268 | C |
---|
269 | CC Method: |
---|
270 | C |
---|
271 | C Declarations: |
---|
272 | C |
---|
273 | |
---|
274 | C |
---|
275 | C Arguments |
---|
276 | INTEGER :: np,nc |
---|
277 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
278 | C |
---|
279 | C Local scalars |
---|
280 | INTEGER :: i,coeffraf,locind_parent_left |
---|
281 | REAL :: ypos,globind_parent_left,globind_parent_right |
---|
282 | REAL :: invds, invds2 |
---|
283 | REAL :: ypos2,diff |
---|
284 | C |
---|
285 | C |
---|
286 | |
---|
287 | coeffraf = nint(ds_parent/ds_child) |
---|
288 | C |
---|
289 | if (coeffraf == 1) then |
---|
290 | C |
---|
291 | return |
---|
292 | C |
---|
293 | endif |
---|
294 | C |
---|
295 | ypos = s_child |
---|
296 | |
---|
297 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
298 | |
---|
299 | globind_parent_left = s_parent |
---|
300 | & + (locind_parent_left - 1)*ds_parent |
---|
301 | |
---|
302 | globind_parent_right = globind_parent_left + ds_parent |
---|
303 | |
---|
304 | C |
---|
305 | invds = 1./ds_parent |
---|
306 | invds2 = ds_child/ds_parent |
---|
307 | |
---|
308 | ypos2 = ypos*invds |
---|
309 | globind_parent_right=globind_parent_right*invds |
---|
310 | |
---|
311 | if (.not.allocated(indparent)) then |
---|
312 | allocate(indparent(nc,1),coeffparent(nc,1)) |
---|
313 | else |
---|
314 | if (size(indparent)<nc) then |
---|
315 | deallocate(indparent,coeffparent) |
---|
316 | allocate(indparent(nc,1),coeffparent(nc,1)) |
---|
317 | endif |
---|
318 | endif |
---|
319 | |
---|
320 | do i = 1,nc-1 |
---|
321 | C |
---|
322 | if (ypos2 > globind_parent_right) then |
---|
323 | locind_parent_left = locind_parent_left + 1 |
---|
324 | globind_parent_right = globind_parent_right + 1. |
---|
325 | endif |
---|
326 | |
---|
327 | diff=(globind_parent_right - ypos2) |
---|
328 | indparent(i,1) = locind_parent_left |
---|
329 | coeffparent(i,1) = diff |
---|
330 | ypos2 = ypos2 + invds2 |
---|
331 | C |
---|
332 | enddo |
---|
333 | C |
---|
334 | ypos = s_child + (nc-1)*ds_child |
---|
335 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
336 | |
---|
337 | if (locind_parent_left == np) then |
---|
338 | indparent(nc,1) = locind_parent_left-1 |
---|
339 | coeffparent(nc,1) = 0. |
---|
340 | else |
---|
341 | C |
---|
342 | globind_parent_left = s_parent |
---|
343 | & + (locind_parent_left - 1)*ds_parent |
---|
344 | C |
---|
345 | indparent(nc,1) = locind_parent_left |
---|
346 | |
---|
347 | coeffparent(nc,1) = (globind_parent_left + ds_parent - ypos) |
---|
348 | & * invds |
---|
349 | endif |
---|
350 | C |
---|
351 | Return |
---|
352 | C |
---|
353 | C |
---|
354 | End Subroutine Linear1dprecompute |
---|
355 | |
---|
356 | Subroutine Linear1daftercompute(x,y,np,nc, |
---|
357 | & s_parent,s_child,ds_parent,ds_child,dir) |
---|
358 | C |
---|
359 | CCC Description: |
---|
360 | CCC Subroutine to do a linear 1D interpolation on a child grid (vector y) from |
---|
361 | CCC its parent grid (vector x) using precomputed coefficient and index. |
---|
362 | C |
---|
363 | CC Method: |
---|
364 | C |
---|
365 | C Declarations: |
---|
366 | C |
---|
367 | |
---|
368 | C |
---|
369 | C Arguments |
---|
370 | INTEGER :: np,nc,dir |
---|
371 | REAL,INTENT(IN), DIMENSION(np) :: x |
---|
372 | REAL,INTENT(OUT), DIMENSION(nc) :: y |
---|
373 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
374 | C |
---|
375 | C Local scalars |
---|
376 | INTEGER :: i,coeffraf,locind_parent_left |
---|
377 | REAL :: ypos,globind_parent_left,globind_parent_right |
---|
378 | REAL :: invds, invds2 |
---|
379 | REAL :: ypos2,diff |
---|
380 | C |
---|
381 | C |
---|
382 | |
---|
383 | !CDIR ALTCODE |
---|
384 | !CDIR NODEP |
---|
385 | do i = 1,nc |
---|
386 | C |
---|
387 | y(i)=coeffparent(i,dir)*x(MAX(indparent(i,dir),1))+ |
---|
388 | & (1.-coeffparent(i,dir))*x(indparent(i,dir)+1) |
---|
389 | C |
---|
390 | enddo |
---|
391 | C |
---|
392 | Return |
---|
393 | C |
---|
394 | C |
---|
395 | End Subroutine Linear1daftercompute |
---|
396 | |
---|
397 | C |
---|
398 | C |
---|
399 | C |
---|
400 | C ************************************************************************** |
---|
401 | CCC Subroutine Lagrange1d |
---|
402 | C ************************************************************************** |
---|
403 | C |
---|
404 | Subroutine Lagrange1d(x,y,np,nc, |
---|
405 | & s_parent,s_child,ds_parent,ds_child) |
---|
406 | C |
---|
407 | CCC Description: |
---|
408 | CCC Subroutine to do a lagrange 1D interpolation on a child grid (vector y) |
---|
409 | CCC from its parent grid (vector x). |
---|
410 | C |
---|
411 | CC Method: |
---|
412 | C |
---|
413 | C Declarations: |
---|
414 | C |
---|
415 | |
---|
416 | C |
---|
417 | C Arguments |
---|
418 | INTEGER :: np,nc |
---|
419 | REAL,INTENT(IN), DIMENSION(np) :: x |
---|
420 | REAL,INTENT(OUT), DIMENSION(nc) :: y |
---|
421 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
422 | C |
---|
423 | C Local scalars |
---|
424 | INTEGER :: i,coeffraf,locind_parent_left |
---|
425 | REAL :: ypos,globind_parent_left |
---|
426 | REAL :: X1,X2,X3 |
---|
427 | real :: deltax,invdsparent |
---|
428 | real t1,t2,t3,t4,t5,t6,t7,t8 |
---|
429 | C |
---|
430 | C |
---|
431 | if (np <= 2) then |
---|
432 | C |
---|
433 | Call Linear1D(x,y,np,nc, |
---|
434 | & s_parent,s_child,ds_parent,ds_child) |
---|
435 | C |
---|
436 | Return |
---|
437 | C |
---|
438 | endif |
---|
439 | C |
---|
440 | coeffraf = nint(ds_parent/ds_child) |
---|
441 | C |
---|
442 | if (coeffraf == 1) then |
---|
443 | C |
---|
444 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
445 | C |
---|
446 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
447 | C |
---|
448 | return |
---|
449 | C |
---|
450 | endif |
---|
451 | |
---|
452 | invdsparent=1./ds_parent |
---|
453 | C |
---|
454 | ypos = s_child |
---|
455 | C |
---|
456 | do i = 1,nc |
---|
457 | C |
---|
458 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
459 | C |
---|
460 | |
---|
461 | globind_parent_left = s_parent |
---|
462 | & + (locind_parent_left - 1)*ds_parent |
---|
463 | |
---|
464 | deltax = invdsparent*(ypos-globind_parent_left) |
---|
465 | ypos = ypos + ds_child |
---|
466 | if (abs(deltax).LE.0.0001) then |
---|
467 | y(i)=x(locind_parent_left) |
---|
468 | |
---|
469 | cycle |
---|
470 | endif |
---|
471 | C |
---|
472 | C |
---|
473 | t2 = deltax - 2. |
---|
474 | t3 = deltax - 1. |
---|
475 | t4 = deltax + 1. |
---|
476 | |
---|
477 | t5 = -(1./6.)*deltax*t2*t3 |
---|
478 | t6 = 0.5*t2*t3*t4 |
---|
479 | t7 = -0.5*deltax*t2*t4 |
---|
480 | t8 = (1./6.)*deltax*t3*t4 |
---|
481 | |
---|
482 | y(i)=t5*x(locind_parent_left-1)+t6*x(locind_parent_left) |
---|
483 | & +t7*x(locind_parent_left+1)+t8*x(locind_parent_left+2) |
---|
484 | C |
---|
485 | C |
---|
486 | enddo |
---|
487 | C |
---|
488 | return |
---|
489 | C |
---|
490 | C |
---|
491 | End Subroutine Lagrange1d |
---|
492 | C |
---|
493 | C |
---|
494 | C ************************************************************************** |
---|
495 | CCC Subroutine Constant1d |
---|
496 | C ************************************************************************** |
---|
497 | C |
---|
498 | Subroutine constant1d(x,y,np,nc, |
---|
499 | & s_parent,s_child,ds_parent,ds_child) |
---|
500 | C |
---|
501 | CCC Description: |
---|
502 | CCC Subroutine to do a linear 1D interpolation on a child grid (vector y) from |
---|
503 | CCC its parent grid (vector x). |
---|
504 | C |
---|
505 | CC Method: |
---|
506 | C |
---|
507 | C Declarations: |
---|
508 | C |
---|
509 | |
---|
510 | C |
---|
511 | C Arguments |
---|
512 | INTEGER :: np,nc |
---|
513 | REAL,INTENT(IN), DIMENSION(np) :: x |
---|
514 | REAL,INTENT(OUT), DIMENSION(nc) :: y |
---|
515 | REAL :: s_parent,s_child,ds_parent,ds_child |
---|
516 | C |
---|
517 | C Local scalars |
---|
518 | INTEGER :: i,coeffraf,locind_parent |
---|
519 | REAL :: ypos |
---|
520 | C |
---|
521 | C |
---|
522 | |
---|
523 | coeffraf = nint(ds_parent/ds_child) |
---|
524 | C |
---|
525 | if (coeffraf == 1) then |
---|
526 | C |
---|
527 | locind_parent = 1 + nint((s_child - s_parent)/ds_parent) |
---|
528 | C |
---|
529 | y(1:nc) = x(locind_parent:locind_parent+nc-1) |
---|
530 | C |
---|
531 | return |
---|
532 | C |
---|
533 | endif |
---|
534 | C |
---|
535 | ypos = s_child |
---|
536 | C |
---|
537 | do i = 1,nc |
---|
538 | C |
---|
539 | locind_parent = 1 + nint((ypos - s_parent)/ds_parent) |
---|
540 | C |
---|
541 | y(i) = x(locind_parent) |
---|
542 | C |
---|
543 | ypos = ypos + ds_child |
---|
544 | C |
---|
545 | enddo |
---|
546 | C |
---|
547 | Return |
---|
548 | C |
---|
549 | C |
---|
550 | End Subroutine constant1d |
---|
551 | C |
---|
552 | C ************************************************************************** |
---|
553 | CCC Subroutine Linear1dconserv |
---|
554 | C ************************************************************************** |
---|
555 | C |
---|
556 | Subroutine Linear1dconserv(x,y,np,nc, |
---|
557 | & s_parent,s_child,ds_parent,ds_child) |
---|
558 | C |
---|
559 | CCC Description: |
---|
560 | CCC Subroutine to do a linear 1D interpolation on a child grid (vector y) from |
---|
561 | CCC its parent grid (vector x). |
---|
562 | C |
---|
563 | CC Method: |
---|
564 | C |
---|
565 | C Declarations: |
---|
566 | C |
---|
567 | Implicit none |
---|
568 | C |
---|
569 | C Arguments |
---|
570 | Integer :: np,nc |
---|
571 | Real, Dimension(np) :: x |
---|
572 | Real, Dimension(nc) :: y |
---|
573 | Real, Dimension(:),Allocatable :: ytemp |
---|
574 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
575 | C |
---|
576 | C Local scalars |
---|
577 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
578 | Real :: ypos |
---|
579 | integer :: i1,i2,ii |
---|
580 | real :: xpmin,xpmax,slope |
---|
581 | INTEGER :: diffmod |
---|
582 | REAL :: xdiffmod |
---|
583 | |
---|
584 | C |
---|
585 | C |
---|
586 | |
---|
587 | coeffraf = nint(ds_parent/ds_child) |
---|
588 | C |
---|
589 | If (coeffraf == 1) Then |
---|
590 | C |
---|
591 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
592 | C |
---|
593 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
594 | C |
---|
595 | return |
---|
596 | C |
---|
597 | End If |
---|
598 | C |
---|
599 | diffmod = 0 |
---|
600 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
601 | |
---|
602 | xdiffmod = real(diffmod)/2. |
---|
603 | |
---|
604 | allocate(ytemp(-2*coeffraf:nc+2*coeffraf)) |
---|
605 | C |
---|
606 | ypos = s_child |
---|
607 | C |
---|
608 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
609 | |
---|
610 | locind_parent_last = 1 + |
---|
611 | & agrif_ceiling((ypos +(nc - 1) *ds_child - s_parent)/ds_parent) |
---|
612 | |
---|
613 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
614 | xpmax = s_parent + (locind_parent_last-1)*ds_parent |
---|
615 | |
---|
616 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
617 | i2 = 1+agrif_int((xpmax-s_child)/ds_child) |
---|
618 | |
---|
619 | i = i1 |
---|
620 | |
---|
621 | if (locind_parent_left == 1) then |
---|
622 | slope= |
---|
623 | & (x(locind_parent_left+1)-x(locind_parent_left))/(coeffraf) |
---|
624 | else |
---|
625 | slope= |
---|
626 | & (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf) |
---|
627 | endif |
---|
628 | |
---|
629 | do ii=i-coeffraf/2+diffmod,i+coeffraf/2 |
---|
630 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
631 | enddo |
---|
632 | |
---|
633 | locind_parent_left = locind_parent_left + 1 |
---|
634 | |
---|
635 | do i=i1 + coeffraf, i2 - coeffraf,coeffraf |
---|
636 | slope= |
---|
637 | & (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf) |
---|
638 | do ii=i-coeffraf/2+diffmod,i+coeffraf/2 |
---|
639 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
640 | enddo |
---|
641 | locind_parent_left = locind_parent_left + 1 |
---|
642 | enddo |
---|
643 | |
---|
644 | i = i2 |
---|
645 | |
---|
646 | if (locind_parent_left == np) then |
---|
647 | slope= |
---|
648 | & (x(locind_parent_left)-x(locind_parent_left-1))/(coeffraf) |
---|
649 | else |
---|
650 | slope= |
---|
651 | & (x(locind_parent_left+1)-x(locind_parent_left-1))/(2.*coeffraf) |
---|
652 | endif |
---|
653 | |
---|
654 | do ii=i-coeffraf/2+diffmod,nc |
---|
655 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
656 | enddo |
---|
657 | C |
---|
658 | y(1:nc)=ytemp(1:nc) |
---|
659 | C |
---|
660 | deallocate(ytemp) |
---|
661 | Return |
---|
662 | C |
---|
663 | End Subroutine Linear1dconserv |
---|
664 | |
---|
665 | C |
---|
666 | C ************************************************************************** |
---|
667 | CCC Subroutine Linear1dconservlim |
---|
668 | C ************************************************************************** |
---|
669 | C |
---|
670 | Subroutine Linear1dconservlim(x,y,np,nc, |
---|
671 | & s_parent,s_child,ds_parent,ds_child) |
---|
672 | C |
---|
673 | CCC Description: |
---|
674 | CCC Subroutine to do a linear 1D interpolation on a child grid (vector y) from |
---|
675 | CCC its parent grid (vector x). |
---|
676 | C |
---|
677 | CC Method: |
---|
678 | C |
---|
679 | C Declarations: |
---|
680 | C |
---|
681 | Implicit none |
---|
682 | C |
---|
683 | C Arguments |
---|
684 | Integer :: np,nc |
---|
685 | Real, Dimension(np) :: x |
---|
686 | Real, Dimension(nc) :: y |
---|
687 | Real, Dimension(:),Allocatable :: ytemp |
---|
688 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
689 | C |
---|
690 | C Local scalars |
---|
691 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
692 | Real :: ypos |
---|
693 | integer :: i1,i2,ii |
---|
694 | real :: xpmin,xpmax,slope |
---|
695 | INTEGER :: diffmod |
---|
696 | real :: xdiffmod |
---|
697 | C |
---|
698 | C |
---|
699 | |
---|
700 | coeffraf = nint(ds_parent/ds_child) |
---|
701 | C |
---|
702 | If (coeffraf == 1) Then |
---|
703 | C |
---|
704 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
705 | C |
---|
706 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
707 | C |
---|
708 | return |
---|
709 | C |
---|
710 | End If |
---|
711 | C |
---|
712 | IF (coeffraf .NE.3) THEN |
---|
713 | print *,'LINEARCONSERVLIM not ready for refinement ratio = ', |
---|
714 | & coeffraf |
---|
715 | stop |
---|
716 | ENDIF |
---|
717 | |
---|
718 | diffmod = 0 |
---|
719 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
720 | |
---|
721 | xdiffmod = real(diffmod)/2. |
---|
722 | |
---|
723 | allocate(ytemp(-2*coeffraf:nc+2*coeffraf)) |
---|
724 | C |
---|
725 | ypos = s_child |
---|
726 | C |
---|
727 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
728 | |
---|
729 | locind_parent_last = 1 + |
---|
730 | & agrif_ceiling((ypos +(nc - 1) *ds_child - s_parent)/ds_parent) |
---|
731 | |
---|
732 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
733 | xpmax = s_parent + (locind_parent_last-1)*ds_parent |
---|
734 | |
---|
735 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
736 | i2 = 1+agrif_int((xpmax-s_child)/ds_child) |
---|
737 | |
---|
738 | i = i1 |
---|
739 | |
---|
740 | if (locind_parent_left == 1) then |
---|
741 | slope=0. |
---|
742 | else |
---|
743 | slope = vanleer(x(locind_parent_left-1:locind_parent_left+1)) |
---|
744 | slope = slope / coeffraf |
---|
745 | endif |
---|
746 | |
---|
747 | do ii=i-coeffraf/2+diffmod,i+coeffraf/2 |
---|
748 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
749 | enddo |
---|
750 | |
---|
751 | locind_parent_left = locind_parent_left + 1 |
---|
752 | |
---|
753 | do i=i1 + coeffraf, i2 - coeffraf,coeffraf |
---|
754 | slope = vanleer(x(locind_parent_left-1:locind_parent_left+1)) |
---|
755 | slope = slope / coeffraf |
---|
756 | |
---|
757 | do ii=i-coeffraf/2+diffmod,i+coeffraf/2 |
---|
758 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
759 | enddo |
---|
760 | locind_parent_left = locind_parent_left + 1 |
---|
761 | enddo |
---|
762 | |
---|
763 | i = i2 |
---|
764 | |
---|
765 | if (locind_parent_left == np) then |
---|
766 | slope=0. |
---|
767 | else |
---|
768 | slope = vanleer(x(locind_parent_left-1:locind_parent_left+1)) |
---|
769 | slope = slope / coeffraf |
---|
770 | endif |
---|
771 | |
---|
772 | do ii=i-coeffraf/2+diffmod,nc |
---|
773 | ytemp(ii) = x(locind_parent_left)+(ii-i-xdiffmod/2.)*slope |
---|
774 | enddo |
---|
775 | C |
---|
776 | y(1:nc)=ytemp(1:nc) |
---|
777 | C |
---|
778 | deallocate(ytemp) |
---|
779 | Return |
---|
780 | C |
---|
781 | End Subroutine Linear1dconservlim |
---|
782 | C |
---|
783 | |
---|
784 | C ************************************************************************** |
---|
785 | CCC Subroutine ppm1d |
---|
786 | C ************************************************************************** |
---|
787 | C |
---|
788 | Subroutine ppm1d(x,y,np,nc, |
---|
789 | & s_parent,s_child,ds_parent,ds_child) |
---|
790 | C |
---|
791 | CCC Description: |
---|
792 | CCC Subroutine to do a 1D interpolation and apply monotonicity constraints |
---|
793 | CCC using piecewise parabolic method |
---|
794 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
795 | CC Method: |
---|
796 | C |
---|
797 | C Declarations: |
---|
798 | C |
---|
799 | Implicit none |
---|
800 | C |
---|
801 | C Arguments |
---|
802 | Integer :: np,nc |
---|
803 | Real, INTENT(IN),Dimension(np) :: x |
---|
804 | Real, INTENT(OUT),Dimension(nc) :: y |
---|
805 | C Real, Dimension(:),Allocatable :: ytemp |
---|
806 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
807 | C |
---|
808 | C Local scalars |
---|
809 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
810 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
811 | Real :: ypos |
---|
812 | integer :: i1,jj |
---|
813 | Real :: xpmin,a |
---|
814 | C |
---|
815 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
816 | Real, Dimension(np) :: xl,delta,a6,slope |
---|
817 | C Real, Dimension(:),Allocatable :: diff,diff2,diff3 |
---|
818 | INTEGER :: diffmod |
---|
819 | REAL :: invcoeffraf |
---|
820 | C |
---|
821 | coeffraf = nint(ds_parent/ds_child) |
---|
822 | C |
---|
823 | If (coeffraf == 1) Then |
---|
824 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
825 | !CDIR ALTCODE |
---|
826 | !CDIR SHORTLOOP |
---|
827 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
828 | return |
---|
829 | End If |
---|
830 | invcoeffraf = ds_child/ds_parent |
---|
831 | C |
---|
832 | |
---|
833 | IF( .NOT. allocated(tabtest4) ) THEN |
---|
834 | Allocate(tabtest4(-2*coeffraf:nc+2*coeffraf)) |
---|
835 | ELSE |
---|
836 | IF (size(tabtest4) .LT. nc+4*coeffraf+1)THEN |
---|
837 | deallocate( tabtest4 ) |
---|
838 | Allocate(tabtest4(-2*coeffraf:nc+2*coeffraf)) |
---|
839 | ENDIF |
---|
840 | ENDIF |
---|
841 | |
---|
842 | ypos = s_child |
---|
843 | C |
---|
844 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
845 | locind_parent_last = 1 + |
---|
846 | & agrif_ceiling((ypos +(nc - 1) |
---|
847 | & *ds_child - s_parent)/ds_parent) |
---|
848 | C |
---|
849 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
850 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
851 | C |
---|
852 | C |
---|
853 | |
---|
854 | !CDIR NOVECTOR |
---|
855 | Do i=1,coeffraf |
---|
856 | tabdiff2(i)=(real(i)-0.5)*invcoeffraf |
---|
857 | EndDo |
---|
858 | |
---|
859 | a = invcoeffraf**2 |
---|
860 | tabdiff3(1) = (1./3.)*a |
---|
861 | a=2.*a |
---|
862 | !CDIR NOVECTOR |
---|
863 | Do i=2,coeffraf |
---|
864 | tabdiff3(i) = tabdiff3(i-1)+(real(i)-1)*a |
---|
865 | EndDo |
---|
866 | C |
---|
867 | if( locind_parent_last+2 <= np ) then |
---|
868 | nmax = locind_parent_last+2 |
---|
869 | else if( locind_parent_last+1 <= np ) then |
---|
870 | nmax = locind_parent_last+1 |
---|
871 | else |
---|
872 | nmax = locind_parent_last |
---|
873 | endif |
---|
874 | C |
---|
875 | if(locind_parent_left-1 >= 1) then |
---|
876 | nmin = locind_parent_left-1 |
---|
877 | else |
---|
878 | nmin = locind_parent_left |
---|
879 | endif |
---|
880 | C |
---|
881 | C |
---|
882 | !CDIR ALTCODE |
---|
883 | !CDIR SHORTLOOP |
---|
884 | Do i = nmin,nmax |
---|
885 | slope(i) = x(i) - x(i-1) |
---|
886 | Enddo |
---|
887 | |
---|
888 | !CDIR ALTCODE |
---|
889 | !CDIR SHORTLOOP |
---|
890 | Do i = nmin+1,nmax-1 |
---|
891 | xl(i)= 0.5*(x(i-1)+x(i)) |
---|
892 | & -0.08333333333333*(slope(i+1)-slope(i-1)) |
---|
893 | Enddo |
---|
894 | C |
---|
895 | C apply parabolic monotonicity |
---|
896 | !CDIR ALTCODE |
---|
897 | !CDIR SHORTLOOP |
---|
898 | Do i = locind_parent_left,locind_parent_last |
---|
899 | delta(i) = xl(i+1) - xl(i) |
---|
900 | a6(i) = 6.*x(i)-3.*(xl(i) +xl(i+1)) |
---|
901 | C |
---|
902 | End do |
---|
903 | C |
---|
904 | diffmod = 0 |
---|
905 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
906 | C |
---|
907 | ipos = i1 |
---|
908 | C |
---|
909 | Do iparent = locind_parent_left,locind_parent_last |
---|
910 | pos=1 |
---|
911 | !CDIR ALTCODE |
---|
912 | !CDIR SHORTLOOP |
---|
913 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
914 | C |
---|
915 | tabtest4(jj) = xl(iparent) |
---|
916 | & + tabdiff2(pos) * (delta(iparent)+a6(iparent)) |
---|
917 | & - tabdiff3(pos) * a6(iparent) |
---|
918 | pos = pos+1 |
---|
919 | End do |
---|
920 | ipos = ipos + coeffraf |
---|
921 | C |
---|
922 | End do |
---|
923 | C |
---|
924 | C |
---|
925 | !CDIR ALTCODE |
---|
926 | !CDIR SHORTLOOP |
---|
927 | y(1:nc)=tabtest4(1:nc) |
---|
928 | |
---|
929 | Return |
---|
930 | End Subroutine ppm1d |
---|
931 | |
---|
932 | |
---|
933 | Subroutine ppm1dprecompute2d(np2,np,nc, |
---|
934 | & s_parent,s_child,ds_parent,ds_child,dir) |
---|
935 | C |
---|
936 | CCC Description: |
---|
937 | CCC Subroutine to compute 2D coefficients and index for a 1D interpolation |
---|
938 | CCC using piecewise parabolic method |
---|
939 | CC Method: |
---|
940 | C |
---|
941 | C Declarations: |
---|
942 | C |
---|
943 | Implicit none |
---|
944 | C |
---|
945 | C Arguments |
---|
946 | Integer :: np2,np,nc,dir |
---|
947 | C Real, Dimension(:),Allocatable :: ytemp |
---|
948 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
949 | C |
---|
950 | C Local scalars |
---|
951 | Integer, Dimension(:,:), Allocatable :: indparent_tmp |
---|
952 | Integer, Dimension(:,:), Allocatable :: indchild_tmp |
---|
953 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
954 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
955 | Real :: ypos |
---|
956 | integer :: i1,jj |
---|
957 | Real :: xpmin,a |
---|
958 | C |
---|
959 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
960 | Real, Dimension(np) :: xl,delta,a6,slope |
---|
961 | INTEGER :: diffmod |
---|
962 | REAL :: invcoeffraf |
---|
963 | C |
---|
964 | coeffraf = nint(ds_parent/ds_child) |
---|
965 | C |
---|
966 | invcoeffraf = ds_child/ds_parent |
---|
967 | C |
---|
968 | |
---|
969 | if (.not.allocated(indparentppm)) then |
---|
970 | allocate( |
---|
971 | & indparentppm(np2*nc,3), |
---|
972 | & indchildppm(np2*nc,3) |
---|
973 | & ) |
---|
974 | else |
---|
975 | if (size(indparentppm,1)<np2*nc) then |
---|
976 | allocate( |
---|
977 | & indparent_tmp(size(indparentppm,1),size(indparentppm,2)), |
---|
978 | & indchild_tmp(size(indparentppm,1),size(indparentppm,2))) |
---|
979 | indparent_tmp = indparentppm |
---|
980 | indchild_tmp = indchildppm |
---|
981 | deallocate(indparentppm,indchildppm) |
---|
982 | allocate( |
---|
983 | &indparentppm(np2*nc,3), |
---|
984 | &indchildppm(np2*nc,3) |
---|
985 | & ) |
---|
986 | indparentppm(1:size(indparent_tmp,1),1:size(indparent_tmp,2)) |
---|
987 | & = indparent_tmp |
---|
988 | indchildppm(1:size(indparent_tmp,1),1:size(indparent_tmp,2)) |
---|
989 | & = indchild_tmp |
---|
990 | deallocate(indparent_tmp,indchild_tmp) |
---|
991 | endif |
---|
992 | endif |
---|
993 | |
---|
994 | if (.not.allocated(indparentppm_1d)) then |
---|
995 | allocate(indparentppm_1d(-2*coeffraf:nc+2*coeffraf), |
---|
996 | & indchildppm_1d(-2*coeffraf:nc+2*coeffraf)) |
---|
997 | else |
---|
998 | if (size(indparentppm_1d)<nc+4*coeffraf+1) then |
---|
999 | deallocate(indparentppm_1d,indchildppm_1d) |
---|
1000 | allocate(indparentppm_1d(-2*coeffraf:nc+2*coeffraf), |
---|
1001 | & indchildppm_1d(-2*coeffraf:nc+2*coeffraf)) |
---|
1002 | endif |
---|
1003 | endif |
---|
1004 | |
---|
1005 | |
---|
1006 | ypos = s_child |
---|
1007 | C |
---|
1008 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
1009 | locind_parent_last = 1 + |
---|
1010 | & agrif_ceiling((ypos +(nc - 1) |
---|
1011 | & *ds_child - s_parent)/ds_parent) |
---|
1012 | C |
---|
1013 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
1014 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
1015 | C |
---|
1016 | C |
---|
1017 | |
---|
1018 | Do i=1,coeffraf |
---|
1019 | tabdiff2(i)=(real(i)-0.5)*invcoeffraf |
---|
1020 | EndDo |
---|
1021 | |
---|
1022 | a = invcoeffraf**2 |
---|
1023 | tabdiff3(1) = (1./3.)*a |
---|
1024 | a=2.*a |
---|
1025 | !CDIR ALTCODE |
---|
1026 | Do i=2,coeffraf |
---|
1027 | tabdiff3(i) = tabdiff3(i-1)+(real(i)-1)*a |
---|
1028 | EndDo |
---|
1029 | |
---|
1030 | !CDIR ALTCODE |
---|
1031 | Do i=1,coeffraf |
---|
1032 | tabppm(1,i,dir) = 0.08333333333333* |
---|
1033 | & (-1.+4*tabdiff2(i)-3*tabdiff3(i)) |
---|
1034 | tabppm(2,i,dir) = 0.08333333333333* |
---|
1035 | & (7.-26.*tabdiff2(i)+18.*tabdiff3(i)) |
---|
1036 | tabppm(3,i,dir) = 0.08333333333333* |
---|
1037 | & (7.+30*tabdiff2(i)-30*tabdiff3(i)) |
---|
1038 | tabppm(4,i,dir) = 0.08333333333333* |
---|
1039 | & (-1.-10.*tabdiff2(i)+18.*tabdiff3(i)) |
---|
1040 | tabppm(5,i,dir) = 0.08333333333333* |
---|
1041 | & (2*tabdiff2(i)-3*tabdiff3(i)) |
---|
1042 | End Do |
---|
1043 | C |
---|
1044 | C |
---|
1045 | diffmod = 0 |
---|
1046 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
1047 | C |
---|
1048 | ipos = i1 |
---|
1049 | C |
---|
1050 | Do iparent = locind_parent_left,locind_parent_last |
---|
1051 | pos=1 |
---|
1052 | !CDIR ALTCODE |
---|
1053 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
1054 | indparentppm_1d(jj) = iparent-2 |
---|
1055 | indchildppm_1d(jj) = pos |
---|
1056 | pos = pos+1 |
---|
1057 | End do |
---|
1058 | ipos = ipos + coeffraf |
---|
1059 | C |
---|
1060 | End do |
---|
1061 | |
---|
1062 | Do i=1,np2 |
---|
1063 | |
---|
1064 | indparentppm(1+(i-1)*nc:i*nc,dir) = |
---|
1065 | & indparentppm_1d(1:nc) + (i-1)*np |
---|
1066 | indchildppm (1+(i-1)*nc:i*nc,dir) = |
---|
1067 | & indchildppm_1d (1:nc) |
---|
1068 | |
---|
1069 | End do |
---|
1070 | |
---|
1071 | Return |
---|
1072 | End Subroutine ppm1dprecompute2d |
---|
1073 | |
---|
1074 | |
---|
1075 | Subroutine ppm1dprecompute(np,nc, |
---|
1076 | & s_parent,s_child,ds_parent,ds_child) |
---|
1077 | C |
---|
1078 | CCC Description: |
---|
1079 | CCC Subroutine to compute coefficient and index for a 1D interpolation |
---|
1080 | CCC using piecewise parabolic method |
---|
1081 | CC Method: |
---|
1082 | C |
---|
1083 | C Declarations: |
---|
1084 | C |
---|
1085 | Implicit none |
---|
1086 | C |
---|
1087 | C Arguments |
---|
1088 | Integer :: np,nc |
---|
1089 | C Real, Dimension(:),Allocatable :: ytemp |
---|
1090 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
1091 | C |
---|
1092 | C Local scalars |
---|
1093 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
1094 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
1095 | Real :: ypos |
---|
1096 | integer :: i1,jj |
---|
1097 | Real :: xpmin,a |
---|
1098 | C |
---|
1099 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
1100 | Real, Dimension(np) :: xl,delta,a6,slope |
---|
1101 | C Real, Dimension(:),Allocatable :: diff,diff2,diff3 |
---|
1102 | INTEGER :: diffmod |
---|
1103 | REAL :: invcoeffraf |
---|
1104 | C |
---|
1105 | coeffraf = nint(ds_parent/ds_child) |
---|
1106 | C |
---|
1107 | If (coeffraf == 1) Then |
---|
1108 | return |
---|
1109 | End If |
---|
1110 | invcoeffraf = ds_child/ds_parent |
---|
1111 | C |
---|
1112 | |
---|
1113 | if (.not.allocated(indparentppm)) then |
---|
1114 | allocate(indparentppm(-2*coeffraf:nc+2*coeffraf,1), |
---|
1115 | & indchildppm(-2*coeffraf:nc+2*coeffraf,1)) |
---|
1116 | else |
---|
1117 | if (size(indparentppm,1)<nc+4*coeffraf+1) then |
---|
1118 | deallocate(indparentppm,indchildppm) |
---|
1119 | allocate(indparentppm(-2*coeffraf:nc+2*coeffraf,1), |
---|
1120 | & indchildppm(-2*coeffraf:nc+2*coeffraf,1)) |
---|
1121 | endif |
---|
1122 | endif |
---|
1123 | |
---|
1124 | ypos = s_child |
---|
1125 | C |
---|
1126 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
1127 | locind_parent_last = 1 + |
---|
1128 | & agrif_ceiling((ypos +(nc - 1) |
---|
1129 | & *ds_child - s_parent)/ds_parent) |
---|
1130 | C |
---|
1131 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
1132 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
1133 | C |
---|
1134 | C |
---|
1135 | |
---|
1136 | Do i=1,coeffraf |
---|
1137 | tabdiff2(i)=(real(i)-0.5)*invcoeffraf |
---|
1138 | EndDo |
---|
1139 | |
---|
1140 | a = invcoeffraf**2 |
---|
1141 | tabdiff3(1) = (1./3.)*a |
---|
1142 | a=2.*a |
---|
1143 | !CDIR ALTCODE |
---|
1144 | !!!CDIR SHORTLOOP |
---|
1145 | Do i=2,coeffraf |
---|
1146 | tabdiff3(i) = tabdiff3(i-1)+(real(i)-1)*a |
---|
1147 | EndDo |
---|
1148 | |
---|
1149 | !CDIR ALTCODE |
---|
1150 | !!!CDIR SHORTLOOP |
---|
1151 | Do i=1,coeffraf |
---|
1152 | tabppm(1,i,1) = 0.08333333333333*(-1.+4*tabdiff2(i)-3*tabdiff3(i)) |
---|
1153 | tabppm(2,i,1) = 0.08333333333333* |
---|
1154 | & (7.-26.*tabdiff2(i)+18.*tabdiff3(i)) |
---|
1155 | tabppm(3,i,1) =0.08333333333333*(7.+30*tabdiff2(i)-30*tabdiff3(i)) |
---|
1156 | tabppm(4,i,1) = 0.08333333333333* |
---|
1157 | & (-1.-10.*tabdiff2(i)+18.*tabdiff3(i)) |
---|
1158 | tabppm(5,i,1) = 0.08333333333333*(2*tabdiff2(i)-3*tabdiff3(i)) |
---|
1159 | End Do |
---|
1160 | C |
---|
1161 | C |
---|
1162 | diffmod = 0 |
---|
1163 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
1164 | C |
---|
1165 | ipos = i1 |
---|
1166 | C |
---|
1167 | Do iparent = locind_parent_left,locind_parent_last |
---|
1168 | pos=1 |
---|
1169 | !CDIR ALTCODE |
---|
1170 | !CDIR SHORTLOOP |
---|
1171 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
1172 | indparentppm(jj,1) = iparent-2 |
---|
1173 | indchildppm(jj,1) = pos |
---|
1174 | pos = pos+1 |
---|
1175 | End do |
---|
1176 | ipos = ipos + coeffraf |
---|
1177 | C |
---|
1178 | End do |
---|
1179 | |
---|
1180 | Return |
---|
1181 | End Subroutine ppm1dprecompute |
---|
1182 | |
---|
1183 | Subroutine ppm1daftercompute(x,y,np,nc, |
---|
1184 | & s_parent,s_child,ds_parent,ds_child,dir) |
---|
1185 | C |
---|
1186 | CCC Description: |
---|
1187 | CCC Subroutine to do a 1D interpolation and apply monotonicity constraints |
---|
1188 | CCC using piecewise parabolic method |
---|
1189 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
1190 | CC Method: Use precomputed coefficient and index |
---|
1191 | C |
---|
1192 | C Declarations: |
---|
1193 | C |
---|
1194 | Implicit none |
---|
1195 | C |
---|
1196 | C Arguments |
---|
1197 | Integer :: np,nc,dir |
---|
1198 | Real, INTENT(IN),Dimension(np) :: x |
---|
1199 | Real, INTENT(OUT),Dimension(nc) :: y |
---|
1200 | C Real, Dimension(:),Allocatable :: ytemp |
---|
1201 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
1202 | C |
---|
1203 | C Local scalars |
---|
1204 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
1205 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
1206 | Real :: ypos |
---|
1207 | integer :: i1,jj |
---|
1208 | Real :: xpmin,a |
---|
1209 | C |
---|
1210 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
1211 | Real, Dimension(np) :: xl,delta,a6,slope |
---|
1212 | INTEGER :: diffmod |
---|
1213 | C |
---|
1214 | C |
---|
1215 | do i=1,nc |
---|
1216 | y(i)=tabppm(1,indchildppm(i,dir),dir)*x(indparentppm(i,dir))+ |
---|
1217 | & tabppm(2,indchildppm(i,dir),dir)*x(indparentppm(i,dir)+1)+ |
---|
1218 | & tabppm(3,indchildppm(i,dir),dir)*x(indparentppm(i,dir)+2)+ |
---|
1219 | & tabppm(4,indchildppm(i,dir),dir)*x(indparentppm(i,dir)+3)+ |
---|
1220 | & tabppm(5,indchildppm(i,dir),dir)*x(indparentppm(i,dir)+4) |
---|
1221 | enddo |
---|
1222 | |
---|
1223 | Return |
---|
1224 | End Subroutine ppm1daftercompute |
---|
1225 | |
---|
1226 | C ************************************************************************** |
---|
1227 | CCC Subroutine weno1d |
---|
1228 | C ************************************************************************** |
---|
1229 | C |
---|
1230 | Subroutine weno1dnew(x,y,np,nc, |
---|
1231 | & s_parent,s_child,ds_parent,ds_child) |
---|
1232 | C |
---|
1233 | CCC Description: |
---|
1234 | CCC Subroutine to do a 1D interpolation and apply monotonicity constraints |
---|
1235 | CCC using piecewise parabolic method |
---|
1236 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
1237 | CC Method: |
---|
1238 | C |
---|
1239 | C Declarations: |
---|
1240 | C |
---|
1241 | Implicit none |
---|
1242 | C |
---|
1243 | C Arguments |
---|
1244 | Integer :: np,nc |
---|
1245 | Real, Dimension(np) :: x |
---|
1246 | Real, Dimension(nc) :: y |
---|
1247 | Real, Dimension(:),Allocatable :: ytemp |
---|
1248 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
1249 | C |
---|
1250 | C Local scalars |
---|
1251 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
1252 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
1253 | Real :: ypos |
---|
1254 | integer :: i1,jj |
---|
1255 | Real :: xpmin,cavg,a,b |
---|
1256 | C |
---|
1257 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
1258 | Real, Dimension(np) :: xr,xl,delta,a6,slope,slope2,smooth |
---|
1259 | Real, Dimension(:),Allocatable :: diff,diff2,diff3 |
---|
1260 | INTEGER :: diffmod |
---|
1261 | REAL :: invcoeffraf |
---|
1262 | integer :: s,l,k |
---|
1263 | integer :: etan, etap |
---|
1264 | real :: delta0, delta1, delta2 |
---|
1265 | real :: epsilon |
---|
1266 | parameter (epsilon = 1.D-8) |
---|
1267 | real, dimension(:,:), allocatable :: ak, ck |
---|
1268 | C |
---|
1269 | coeffraf = nint(ds_parent/ds_child) |
---|
1270 | C |
---|
1271 | If (coeffraf == 1) Then |
---|
1272 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
1273 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
1274 | return |
---|
1275 | End If |
---|
1276 | invcoeffraf = ds_child/ds_parent |
---|
1277 | Allocate(ak(0:1,coeffraf)) |
---|
1278 | Allocate(ck(0:1,coeffraf)) |
---|
1279 | |
---|
1280 | C |
---|
1281 | Allocate(ytemp(-2*coeffraf:nc+2*coeffraf)) |
---|
1282 | ypos = s_child |
---|
1283 | C |
---|
1284 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
1285 | locind_parent_last = 1 + |
---|
1286 | & agrif_ceiling((ypos +(nc - 1) |
---|
1287 | & *ds_child - s_parent)/ds_parent) |
---|
1288 | C |
---|
1289 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
1290 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
1291 | C |
---|
1292 | Allocate( diff(coeffraf),diff2(coeffraf),diff3(coeffraf) ) |
---|
1293 | C |
---|
1294 | diff(1)=0.5*invcoeffraf |
---|
1295 | do i=2,coeffraf |
---|
1296 | diff(i) = diff(i-1)+invcoeffraf |
---|
1297 | enddo |
---|
1298 | |
---|
1299 | ak = 0. |
---|
1300 | ck = 0. |
---|
1301 | |
---|
1302 | do i=1,coeffraf |
---|
1303 | do k=0,1 |
---|
1304 | do s=0,2 |
---|
1305 | do l=0,2 |
---|
1306 | if (l /= s) then |
---|
1307 | ak(k,i) = ak(k,i)+(diff(i)-(k-l+1.)) |
---|
1308 | endif |
---|
1309 | enddo |
---|
1310 | enddo |
---|
1311 | enddo |
---|
1312 | |
---|
1313 | etap = 0 |
---|
1314 | etan = 0 |
---|
1315 | do k=0,1 |
---|
1316 | if (ak(k,i) > 0) then |
---|
1317 | etap = etap+1 |
---|
1318 | else if (ak(k,i) < 0) then |
---|
1319 | etan = etan + 1 |
---|
1320 | endif |
---|
1321 | enddo |
---|
1322 | |
---|
1323 | do k=0,1 |
---|
1324 | if (ak(k,i) == 0) then |
---|
1325 | Ck(k,i) = 1. |
---|
1326 | else if (ak(k,i) > 0) then |
---|
1327 | Ck(k,i) = 1./(etap * ak(k,i)) |
---|
1328 | else |
---|
1329 | Ck(k,i) = -1./(etan * ak(k,i)) |
---|
1330 | endif |
---|
1331 | enddo |
---|
1332 | enddo |
---|
1333 | |
---|
1334 | C |
---|
1335 | a = 0. |
---|
1336 | b = invcoeffraf |
---|
1337 | Do i=1,coeffraf |
---|
1338 | diff2(i) = 0.5*(b*b - a*a) |
---|
1339 | diff3(i) = (1./3.)*(b*b*b - a*a*a) |
---|
1340 | a = a + invcoeffraf |
---|
1341 | b = b + invcoeffraf |
---|
1342 | End do |
---|
1343 | C |
---|
1344 | if( locind_parent_last+2 <= np ) then |
---|
1345 | nmax = locind_parent_last+2 |
---|
1346 | elseif( locind_parent_last+1 <= np ) then |
---|
1347 | nmax = locind_parent_last+1 |
---|
1348 | else |
---|
1349 | nmax = locind_parent_last |
---|
1350 | endif |
---|
1351 | C |
---|
1352 | if(locind_parent_left-2 >= 1) then |
---|
1353 | nmin = locind_parent_left-2 |
---|
1354 | elseif(locind_parent_left-1 >= 1) then |
---|
1355 | nmin = locind_parent_left-1 |
---|
1356 | else |
---|
1357 | nmin = locind_parent_left |
---|
1358 | endif |
---|
1359 | C |
---|
1360 | Do i = nmin+1,nmax |
---|
1361 | slope(i) = (x(i) - x(i-1)) |
---|
1362 | Enddo |
---|
1363 | DO i=nmin+2,nmax |
---|
1364 | smooth(i) = 0.5*(slope(i)**2+slope(i-1)**2) |
---|
1365 | & +(slope(i)-slope(i-1))**2 |
---|
1366 | enddo |
---|
1367 | C |
---|
1368 | diffmod = 0 |
---|
1369 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
1370 | C |
---|
1371 | ipos = i1 |
---|
1372 | C |
---|
1373 | Do iparent = locind_parent_left,locind_parent_last |
---|
1374 | pos=1 |
---|
1375 | |
---|
1376 | delta0=1./(epsilon+smooth(iparent ))**3 |
---|
1377 | delta1=1./(epsilon+smooth(iparent+1))**3 |
---|
1378 | delta2=1./(epsilon+smooth(iparent+2))**3 |
---|
1379 | |
---|
1380 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
1381 | C |
---|
1382 | pos = pos+1 |
---|
1383 | End do |
---|
1384 | ipos = ipos + coeffraf |
---|
1385 | C |
---|
1386 | End do |
---|
1387 | C |
---|
1388 | C |
---|
1389 | y(1:nc)=ytemp(1:nc) |
---|
1390 | deallocate(ytemp) |
---|
1391 | deallocate(diff, diff2, diff3) |
---|
1392 | |
---|
1393 | deallocate(ak,ck) |
---|
1394 | |
---|
1395 | Return |
---|
1396 | End Subroutine weno1dnew |
---|
1397 | |
---|
1398 | C ************************************************************************** |
---|
1399 | CCC Subroutine weno1d |
---|
1400 | C ************************************************************************** |
---|
1401 | C |
---|
1402 | Subroutine weno1d(x,y,np,nc, |
---|
1403 | & s_parent,s_child,ds_parent,ds_child) |
---|
1404 | C |
---|
1405 | CCC Description: |
---|
1406 | CCC Subroutine to do a 1D interpolation and apply monotonicity constraints |
---|
1407 | CCC using piecewise parabolic method |
---|
1408 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
1409 | CC Method: |
---|
1410 | C |
---|
1411 | C Declarations: |
---|
1412 | C |
---|
1413 | Implicit none |
---|
1414 | C |
---|
1415 | C Arguments |
---|
1416 | Integer :: np,nc |
---|
1417 | Real, Dimension(np) :: x |
---|
1418 | Real, Dimension(nc) :: y |
---|
1419 | Real, Dimension(:),Allocatable :: ytemp |
---|
1420 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
1421 | C |
---|
1422 | C Local scalars |
---|
1423 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
1424 | Integer :: iparent,ipos,pos,nmin,nmax |
---|
1425 | Real :: ypos |
---|
1426 | integer :: i1,jj |
---|
1427 | Real :: xpmin,cavg,a,b |
---|
1428 | C |
---|
1429 | Real :: xrmin,xrmax,am3,s2,s1 |
---|
1430 | Real, Dimension(np) :: xr,xl,delta,a6,slope,slope2 |
---|
1431 | Real, Dimension(:),Allocatable :: diff,diff2,diff3 |
---|
1432 | INTEGER :: diffmod |
---|
1433 | REAL :: invcoeffraf |
---|
1434 | integer :: s,l,k |
---|
1435 | integer :: etan, etap |
---|
1436 | real :: delta0, delta1,sumdelta |
---|
1437 | real :: epsilon |
---|
1438 | parameter (epsilon = 1.D-8) |
---|
1439 | C |
---|
1440 | coeffraf = nint(ds_parent/ds_child) |
---|
1441 | C |
---|
1442 | If (coeffraf == 1) Then |
---|
1443 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
1444 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
1445 | return |
---|
1446 | End If |
---|
1447 | invcoeffraf = ds_child/ds_parent |
---|
1448 | |
---|
1449 | C |
---|
1450 | Allocate(ytemp(-2*coeffraf:nc+2*coeffraf)) |
---|
1451 | ypos = s_child |
---|
1452 | C |
---|
1453 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
1454 | locind_parent_last = 1 + |
---|
1455 | & agrif_ceiling((ypos +(nc - 1) |
---|
1456 | & *ds_child - s_parent)/ds_parent) |
---|
1457 | C |
---|
1458 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
1459 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
1460 | C |
---|
1461 | Allocate( diff(coeffraf)) |
---|
1462 | C |
---|
1463 | diff(1)=0.5*invcoeffraf |
---|
1464 | do i=2,coeffraf |
---|
1465 | diff(i) = diff(i-1)+invcoeffraf |
---|
1466 | enddo |
---|
1467 | C |
---|
1468 | if( locind_parent_last+2 <= np ) then |
---|
1469 | nmax = locind_parent_last+2 |
---|
1470 | else if( locind_parent_last+1 <= np ) then |
---|
1471 | nmax = locind_parent_last+1 |
---|
1472 | else |
---|
1473 | nmax = locind_parent_last |
---|
1474 | endif |
---|
1475 | C |
---|
1476 | if(locind_parent_left-1 >= 1) then |
---|
1477 | nmin = locind_parent_left-1 |
---|
1478 | else |
---|
1479 | nmin = locind_parent_left |
---|
1480 | endif |
---|
1481 | C |
---|
1482 | Do i = nmin+1,nmax |
---|
1483 | slope(i) = (x(i) - x(i-1)) |
---|
1484 | Enddo |
---|
1485 | C |
---|
1486 | diffmod = 0 |
---|
1487 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
1488 | C |
---|
1489 | ipos = i1 |
---|
1490 | C |
---|
1491 | Do iparent = locind_parent_left,locind_parent_last |
---|
1492 | pos=1 |
---|
1493 | delta0=1./(epsilon+slope(iparent )**2)**2 |
---|
1494 | delta1=1./(epsilon+slope(iparent+1)**2)**2 |
---|
1495 | sumdelta = 1./(delta0+delta1) |
---|
1496 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
1497 | C |
---|
1498 | ytemp(jj) = x(iparent)+(diff(pos)-0.5)*( |
---|
1499 | & delta0*slope(iparent)+ |
---|
1500 | & delta1*slope(iparent+1))*sumdelta |
---|
1501 | pos = pos+1 |
---|
1502 | End do |
---|
1503 | ipos = ipos + coeffraf |
---|
1504 | C |
---|
1505 | End do |
---|
1506 | C |
---|
1507 | C |
---|
1508 | y(1:nc)=ytemp(1:nc) |
---|
1509 | deallocate(ytemp) |
---|
1510 | deallocate(diff) |
---|
1511 | |
---|
1512 | Return |
---|
1513 | End Subroutine weno1d |
---|
1514 | |
---|
1515 | C |
---|
1516 | C ************************************************************************** |
---|
1517 | CCC Subroutine eno1d |
---|
1518 | C ************************************************************************** |
---|
1519 | C |
---|
1520 | Subroutine eno1d(x,y,np,nc, |
---|
1521 | & s_parent,s_child,ds_parent,ds_child) |
---|
1522 | C |
---|
1523 | CCC Description: |
---|
1524 | CCC ---- p 163-164 Computational gasdynamics ---- |
---|
1525 | CCC Subroutine to do a 1D interpolation |
---|
1526 | CCC using piecewise polynomial ENO reconstruction technique |
---|
1527 | CCC on a child grid (vector y) from its parent grid (vector x). |
---|
1528 | CC Method: |
---|
1529 | C |
---|
1530 | C Declarations: |
---|
1531 | C |
---|
1532 | Implicit none |
---|
1533 | C |
---|
1534 | C Arguments |
---|
1535 | Integer :: np,nc |
---|
1536 | Real, Dimension(np) :: x |
---|
1537 | Real, Dimension(nc) :: y |
---|
1538 | Real, Dimension(:),Allocatable :: ytemp |
---|
1539 | Real :: s_parent,s_child,ds_parent,ds_child |
---|
1540 | C |
---|
1541 | C Local scalars |
---|
1542 | Integer :: i,coeffraf,locind_parent_left,locind_parent_last |
---|
1543 | Integer :: ipos,pos |
---|
1544 | Real :: ypos,xi |
---|
1545 | integer :: i1,jj |
---|
1546 | Real :: xpmin,cavg |
---|
1547 | C |
---|
1548 | Real, Dimension(3,np) :: dd,c |
---|
1549 | Integer :: left |
---|
1550 | C |
---|
1551 | Real, DImension(1:np+1) :: xhalf |
---|
1552 | Real, Dimension(:,:),Allocatable :: Xbar |
---|
1553 | INTEGER :: diffmod |
---|
1554 | C |
---|
1555 | coeffraf = nint(ds_parent/ds_child) |
---|
1556 | C |
---|
1557 | If (coeffraf == 1) Then |
---|
1558 | locind_parent_left = 1 + nint((s_child - s_parent)/ds_parent) |
---|
1559 | y(1:nc) = x(locind_parent_left:locind_parent_left+nc-1) |
---|
1560 | return |
---|
1561 | End If |
---|
1562 | |
---|
1563 | diffmod = 0 |
---|
1564 | IF (mod(coeffraf,2) == 0) diffmod = 1 |
---|
1565 | C |
---|
1566 | Allocate(ytemp(-2*coeffraf:nc+2*coeffraf)) |
---|
1567 | ypos = s_child |
---|
1568 | locind_parent_left = 1 + agrif_int((ypos - s_parent)/ds_parent) |
---|
1569 | locind_parent_last = 1 + |
---|
1570 | & agrif_ceiling((ypos +(nc - 1) *ds_child - |
---|
1571 | & s_parent)/ds_parent) |
---|
1572 | xpmin = s_parent + (locind_parent_left-1)*ds_parent |
---|
1573 | i1 = 1+agrif_int((xpmin-s_child)/ds_child) |
---|
1574 | C |
---|
1575 | xhalf(np+1) = np + 0.5 |
---|
1576 | Do i = 1,np |
---|
1577 | xhalf(i) = i - 0.5 |
---|
1578 | Enddo |
---|
1579 | C |
---|
1580 | C compute divided differences |
---|
1581 | C |
---|
1582 | dd(1,1:np) = x(1:np) |
---|
1583 | dd(2,1:np-1) = 0.5*( dd(1,2:np) - dd(1,1:np-1) ) |
---|
1584 | dd(3,1:np-2) = (1./3.)*( dd(2,2:np-1) - dd(2,1:np-2) ) |
---|
1585 | C |
---|
1586 | Allocate( Xbar( coeffraf,2 ) ) |
---|
1587 | xi = 0.5 |
---|
1588 | Do i = 1,coeffraf |
---|
1589 | Xbar(i,1) = (i-1)*ds_child/ds_parent - xi |
---|
1590 | Xbar(i,2) = i*ds_child/ds_parent - xi |
---|
1591 | Enddo |
---|
1592 | C |
---|
1593 | ipos = i1 |
---|
1594 | C |
---|
1595 | DO i = locind_parent_left,locind_parent_last |
---|
1596 | left = i |
---|
1597 | do jj = 2,3 |
---|
1598 | If(abs(dd(jj,left)) .gt. abs(dd(jj,left-1))) |
---|
1599 | & left = left-1 |
---|
1600 | enddo |
---|
1601 | C |
---|
1602 | C convert to Taylor series form |
---|
1603 | C |
---|
1604 | Call Taylor(i,xhalf(left:left+2),dd(1:3,left),c(1:3,i)) |
---|
1605 | ENDDO |
---|
1606 | C |
---|
1607 | C evaluate the reconstruction on each cell |
---|
1608 | C |
---|
1609 | DO i = locind_parent_left,locind_parent_last |
---|
1610 | C |
---|
1611 | cavg = 0. |
---|
1612 | pos = 1. |
---|
1613 | C |
---|
1614 | Do jj = ipos - coeffraf/2+diffmod,ipos + coeffraf/2 |
---|
1615 | ytemp(jj) =(c(1,i)*(Xbar(pos,2)-Xbar(pos,1)) |
---|
1616 | & +c(2,i)*(Xbar(pos,2)*Xbar(pos,2)- |
---|
1617 | & Xbar(pos,1)*Xbar(pos,1)) |
---|
1618 | & +c(3,i)*(Xbar(pos,2)*Xbar(pos,2)*Xbar(pos,2)- |
---|
1619 | & Xbar(pos,1)*Xbar(pos,1)*Xbar(pos,1))) |
---|
1620 | & *coeffraf |
---|
1621 | cavg = cavg + ytemp(jj) |
---|
1622 | pos = pos+1 |
---|
1623 | Enddo |
---|
1624 | ipos = ipos + coeffraf |
---|
1625 | ENDDO |
---|
1626 | C |
---|
1627 | y(1:nc)=ytemp(1:nc) |
---|
1628 | deallocate(ytemp,Xbar) |
---|
1629 | C |
---|
1630 | Return |
---|
1631 | End Subroutine eno1d |
---|
1632 | C |
---|
1633 | C |
---|
1634 | C ************************************************************************** |
---|
1635 | CCC Subroutine taylor |
---|
1636 | C ************************************************************************** |
---|
1637 | C |
---|
1638 | subroutine taylor(ind,xhalf,dd,c) |
---|
1639 | C |
---|
1640 | Integer :: ind |
---|
1641 | real,dimension(3) :: dd,c |
---|
1642 | real,dimension(0:3,0:3) :: d |
---|
1643 | real,dimension(3) :: xhalf |
---|
1644 | integer ::i,j |
---|
1645 | C |
---|
1646 | C |
---|
1647 | d(0,0:3)=1. |
---|
1648 | do i = 1,3 |
---|
1649 | d(i,0)=(ind-xhalf(i))*d(i-1,0) |
---|
1650 | enddo |
---|
1651 | C |
---|
1652 | do i = 1,3 |
---|
1653 | do j = 1,3-i |
---|
1654 | d(i,j) = d(i,j-1) + (ind-xhalf(i+j))*d(i-1,j) |
---|
1655 | enddo |
---|
1656 | enddo |
---|
1657 | C |
---|
1658 | do j = 1,3 |
---|
1659 | c(j) = 0. |
---|
1660 | do i=0,3-j |
---|
1661 | c(j) = c(j) + d(i,j)*dd(i+j) |
---|
1662 | enddo |
---|
1663 | enddo |
---|
1664 | C |
---|
1665 | end subroutine taylor |
---|
1666 | |
---|
1667 | |
---|
1668 | REAL FUNCTION vanleer(tab) |
---|
1669 | REAL, DIMENSION(3) :: tab |
---|
1670 | real res1 |
---|
1671 | real p1,p2,p3 |
---|
1672 | |
---|
1673 | p1=(tab(3)-tab(1))/2. |
---|
1674 | p2=2.*(tab(2)-tab(1)) |
---|
1675 | p3=2.*(tab(3)-tab(2)) |
---|
1676 | |
---|
1677 | if ((p1>0.).AND.(p2>0.).AND.(p3>0)) then |
---|
1678 | res1=minval((/p1,p2,p3/)) |
---|
1679 | elseif ((p1<0.).AND.(p2<0.).AND.(p3<0)) then |
---|
1680 | res1=maxval((/p1,p2,p3/)) |
---|
1681 | else |
---|
1682 | res1=0. |
---|
1683 | endif |
---|
1684 | |
---|
1685 | vanleer = res1 |
---|
1686 | |
---|
1687 | |
---|
1688 | END FUNCTION vanleer |
---|
1689 | |
---|
1690 | C |
---|
1691 | End Module Agrif_Interpbasic |
---|