17 implicit none ;
private 19 #include <MOM_memory.h> 25 integer :: remapping_scheme = -911
29 logical :: boundary_extrapolation = .true.
31 logical :: check_reconstruction = .false.
33 logical :: check_remapping = .false.
35 logical :: force_bounds_in_subcell = .false.
57 character(len=40) ::
mdl =
"MOM_remapping" 61 "PCM (1st-order accurate)\n"//&
62 "PLM (2nd-order accurate)\n"//&
63 "PPM_H4 (3rd-order accurate)\n"//&
64 "PPM_IH4 (3rd-order accurate)\n"//&
65 "PQM_IH4IH3 (4th-order accurate)\n"//&
66 "PQM_IH6IH5 (5th-order accurate)\n" 72 #define __USE_ROUNDOFF_SAFE_ADJUSTMENTS__ 87 check_reconstruction, check_remapping, force_bounds_in_subcell)
89 character(len=*),
optional,
intent(in) :: remapping_scheme
90 logical,
optional,
intent(in) :: boundary_extrapolation
91 logical,
optional,
intent(in) :: check_reconstruction
92 logical,
optional,
intent(in) :: check_remapping
93 logical,
optional,
intent(in) :: force_bounds_in_subcell
95 if (
present(remapping_scheme))
then 98 if (
present(boundary_extrapolation))
then 99 cs%boundary_extrapolation = boundary_extrapolation
101 if (
present(check_reconstruction))
then 102 cs%check_reconstruction = check_reconstruction
104 if (
present(check_remapping))
then 105 cs%check_remapping = check_remapping
107 if (
present(force_bounds_in_subcell))
then 108 cs%force_bounds_in_subcell = force_bounds_in_subcell
114 integer,
intent(in) :: nz
115 real,
dimension(nz),
intent(in) :: h
116 real,
dimension(nz+1),
intent(inout) :: x
136 integer,
intent(in) :: n1
137 integer,
intent(in) :: n2
138 real,
intent(in) :: sum1
139 real,
intent(in) :: sum2
140 logical :: isPosSumErrSignificant
142 real :: sumErr, allowedErr, eps
144 if (sum1<0.)
call mom_error(fatal,
'isPosSumErrSignificant: sum1<0 is not allowed!')
145 if (sum2<0.)
call mom_error(fatal,
'isPosSumErrSignificant: sum2<0 is not allowed!')
146 sumerr = abs(sum1-sum2)
148 allowederr = eps*0.5*(
real(n1-1)*sum1+
real(n2-1)*sum2)
149 if (sumerr>allowederr)
then 150 write(0,*)
'isPosSumErrSignificant: sum1,sum2=',sum1,sum2
151 write(0,*)
'isPosSumErrSignificant: eps=',eps
152 write(0,*)
'isPosSumErrSignificant: err,n*eps=',sumerr,allowederr
153 write(0,*)
'isPosSumErrSignificant: err/eps,n1,n2,n1+n2=',sumerr/eps,n1,n2,n1+n2
154 ispossumerrsignificant = .true.
156 ispossumerrsignificant = .false.
164 integer,
intent(in) :: n0
165 real,
dimension(n0),
intent(in) :: h0
166 real,
dimension(n0),
intent(in) :: u0
167 integer,
intent(in) :: n1
168 real,
dimension(n1),
intent(in) :: h1
169 real,
dimension(n1),
intent(out) :: u1
172 real,
dimension(n0,2) :: ppoly_r_E
173 real,
dimension(n0,2) :: ppoly_r_S
174 real,
dimension(n0,CS%degree+1) :: ppoly_r_coefficients
176 real :: eps, h0tot, h0err, h1tot, h1err, u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, uh_err
179 ppoly_r_coefficients, ppoly_r_e, ppoly_r_s, imethod )
182 cs%boundary_extrapolation, ppoly_r_coefficients, ppoly_r_e, ppoly_r_s)
186 cs%force_bounds_in_subcell, u1, uh_err )
188 if (cs%check_remapping)
then 193 if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err .and. abs(h1tot-h0tot)<h0err+h1err) &
194 .or. (u1min<u0min .or. u1max>u0max) )
then 195 write(0,*)
'iMethod = ',imethod
196 write(0,*)
'H: h0tot=',h0tot,
'h1tot=',h1tot,
'dh=',h1tot-h0tot,
'h0err=',h0err,
'h1err=',h1err
197 if (abs(h1tot-h0tot)>h0err+h1err)
write(0,*)
'H non-conservation difference=',h1tot-h0tot,
'allowed err=',h0err+h1err,
' <-----!' 198 write(0,*)
'UH: u0tot=',u0tot,
'u1tot=',u1tot,
'duh=',u1tot-u0tot,
'u0err=',u0err,
'u1err=',u1err,
'uh_err=',uh_err
199 if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err)
write(0,*)
'U non-conservation difference=',u1tot-u0tot,
'allowed err=',u0err+u1err+uh_err,
' <-----!' 200 write(0,*)
'U: u0min=',u0min,
'u1min=',u1min
201 if (u1min<u0min)
write(0,*)
'U minimum overshoot=',u1min-u0min,
' <-----!' 202 write(0,*)
'U: u0max=',u0max,
'u1max=',u1max
203 if (u1max>u0max)
write(0,*)
'U maximum overshoot=',u1max-u0max,
' <-----!' 204 write(0,
'(a3,6a24)')
'k',
'h0',
'left edge',
'u0',
'right edge',
'h1',
'u1' 206 if (k<=min(n0,n1))
then 207 write(0,
'(i3,1p6e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2),h1(k),u1(k)
209 write(0,
'(i3,96x,1p2e24.16)') k,h1(k),u1(k)
211 write(0,
'(i3,1p4e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2)
214 write(0,
'(a3,2a24)')
'k',
'u0',
'Polynomial coefficients' 216 write(0,
'(i3,1p6e24.16)') k,u0(k),ppoly_r_coefficients(k,:)
218 call mom_error( fatal,
'MOM_remapping, remapping_core_h: '//&
219 'Remapping result is inconsistent!' )
230 integer,
intent(in) :: n0
231 real,
dimension(n0),
intent(in) :: h0
232 real,
dimension(n0),
intent(in) :: u0
233 integer,
intent(in) :: n1
234 real,
dimension(n1+1),
intent(in) :: dx
235 real,
dimension(n1),
intent(out) :: u1
238 real,
dimension(n0,2) :: ppoly_r_E
239 real,
dimension(n0,2) :: ppoly_r_S
240 real,
dimension(n0,CS%degree+1) :: ppoly_r_coefficients
242 real :: eps, h0tot, h0err, h1tot, h1err
243 real :: u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, uh_err
244 real,
dimension(n1) :: h1
247 ppoly_r_coefficients, ppoly_r_e, ppoly_r_s, imethod )
250 cs%boundary_extrapolation, ppoly_r_coefficients, ppoly_r_e, ppoly_r_s)
255 h1(k) = max( 0., h0(k) + ( dx(k+1) - dx(k) ) )
257 h1(k) = max( 0., dx(k+1) - dx(k) )
261 cs%force_bounds_in_subcell,u1, uh_err )
265 if (cs%check_remapping)
then 270 if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err .and. abs(h1tot-h0tot)<h0err+h1err) &
271 .or. (u1min<u0min .or. u1max>u0max) )
then 272 write(0,*)
'iMethod = ',imethod
273 write(0,*)
'H: h0tot=',h0tot,
'h1tot=',h1tot,
'dh=',h1tot-h0tot,
'h0err=',h0err,
'h1err=',h1err
274 if (abs(h1tot-h0tot)>h0err+h1err)
write(0,*)
'H non-conservation difference=',h1tot-h0tot,
'allowed err=',h0err+h1err,
' <-----!' 275 write(0,*)
'UH: u0tot=',u0tot,
'u1tot=',u1tot,
'duh=',u1tot-u0tot,
'u0err=',u0err,
'u1err=',u1err,
'uh_err=',uh_err
276 if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err)
write(0,*)
'U non-conservation difference=',u1tot-u0tot,
'allowed err=',u0err+u1err+uh_err,
' <-----!' 277 write(0,*)
'U: u0min=',u0min,
'u1min=',u1min
278 if (u1min<u0min)
write(0,*)
'U minimum overshoot=',u1min-u0min,
' <-----!' 279 write(0,*)
'U: u0max=',u0max,
'u1max=',u1max
280 if (u1max>u0max)
write(0,*)
'U maximum overshoot=',u1max-u0max,
' <-----!' 281 write(0,
'(a3,6a24)')
'k',
'h0',
'left edge',
'u0',
'right edge',
'h1',
'u1' 283 if (k<=min(n0,n1))
then 284 write(0,
'(i3,1p6e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2),h1(k),u1(k)
286 write(0,
'(i3,96x,1p2e24.16)') k,h1(k),u1(k)
288 write(0,
'(i3,1p4e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2)
291 write(0,
'(a3,2a24)')
'k',
'u0',
'Polynomial coefficients' 293 write(0,
'(i3,1p6e24.16)') k,u0(k),ppoly_r_coefficients(k,:)
295 call mom_error( fatal,
'MOM_remapping, remapping_core_w: '//&
296 'Remapping result is inconsistent!' )
305 ppoly_r_coefficients, ppoly_r_E, ppoly_r_S, iMethod )
306 integer,
intent(in) :: n0
307 real,
dimension(n0),
intent(in) :: h0
308 real,
dimension(n0),
intent(in) :: u0
309 integer,
intent(in) :: remapping_scheme
310 integer,
intent(in) :: deg
311 logical,
intent(in) :: boundary_extrapolation
312 real,
dimension(n0,deg+1),
intent(out) :: ppoly_r_coefficients
313 real,
dimension(n0,2),
intent(out) :: ppoly_r_E
314 real,
dimension(n0,2),
intent(out) :: ppoly_r_S
315 integer,
intent(out) :: iMethod
317 integer :: local_remapping_scheme
322 ppoly_r_coefficients(:,:) = 0.0
325 local_remapping_scheme = remapping_scheme
329 local_remapping_scheme = min( local_remapping_scheme,
remapping_plm )
333 select case ( local_remapping_scheme )
339 if ( boundary_extrapolation)
then 344 call edge_values_explicit_h4( n0, h0, u0, ppoly_r_e )
346 if ( boundary_extrapolation)
then 353 if ( boundary_extrapolation)
then 361 if ( boundary_extrapolation)
then 369 if ( boundary_extrapolation)
then 374 call mom_error( fatal,
'MOM_remapping, build_reconstructions_1d: '//&
375 'The selected remapping method is invalid' )
382 ppoly_r_coefficients, ppoly_r_E, ppoly_r_S)
383 integer,
intent(in) :: n0
384 real,
dimension(n0),
intent(in) :: h0
385 real,
dimension(n0),
intent(in) :: u0
386 integer,
intent(in) :: deg
387 logical,
intent(in) :: boundary_extrapolation
388 real,
dimension(n0,deg+1),
intent(out) :: ppoly_r_coefficients
389 real,
dimension(n0,2),
intent(out) :: ppoly_r_E
390 real,
dimension(n0,2),
intent(out) :: ppoly_r_S
393 real :: u_l, u_c, u_r
395 logical :: problem_detected
397 problem_detected = .false.
399 u_l = u0(max(1,i0-1))
401 u_r = u0(min(n0,i0+1))
402 if (i0 > 1 .or. .not. boundary_extrapolation)
then 403 u_min = min(u_l, u_c)
404 u_max = max(u_l, u_c)
405 if (ppoly_r_e(i0,1) < u_min)
then 406 write(0,
'(a,i4,5(x,a,1pe24.16))')
'Left edge undershoot at',i0,
'u(i0-1)=',u_l,
'u(i0)=',u_c, &
407 'edge=',ppoly_r_e(i0,1),
'err=',ppoly_r_e(i0,1)-u_min
408 problem_detected = .true.
410 if (ppoly_r_e(i0,1) > u_max)
then 411 write(0,
'(a,i4,5(x,a,1pe24.16))')
'Left edge overshoot at',i0,
'u(i0-1)=',u_l,
'u(i0)=',u_c, &
412 'edge=',ppoly_r_e(i0,1),
'err=',ppoly_r_e(i0,1)-u_max
413 problem_detected = .true.
416 if (i0 < n0 .or. .not. boundary_extrapolation)
then 417 u_min = min(u_c, u_r)
418 u_max = max(u_c, u_r)
419 if (ppoly_r_e(i0,2) < u_min)
then 420 write(0,
'(a,i4,5(x,a,1pe24.16))')
'Right edge undershoot at',i0,
'u(i0)=',u_c,
'u(i0+1)=',u_r, &
421 'edge=',ppoly_r_e(i0,2),
'err=',ppoly_r_e(i0,2)-u_min
422 problem_detected = .true.
424 if (ppoly_r_e(i0,2) > u_max)
then 425 write(0,
'(a,i4,5(x,a,1pe24.16))')
'Right edge overshoot at',i0,
'u(i0)=',u_c,
'u(i0+1)=',u_r, &
426 'edge=',ppoly_r_e(i0,2),
'err=',ppoly_r_e(i0,2)-u_max
427 problem_detected = .true.
431 if ( (u_c-u_l)*(ppoly_r_e(i0,1)-ppoly_r_e(i0-1,2)) < 0.)
then 432 write(0,
'(a,i4,5(x,a,1pe24.16))')
'Non-monotonic edges at',i0,
'u(i0-1)=',u_l,
'u(i0)=',u_c, &
433 'right edge=',ppoly_r_e(i0-1,2),
'left edge=',ppoly_r_e(i0,1)
434 write(0,
'(5(a,1pe24.16,x))')
'u(i0)-u(i0-1)',u_c-u_l,
'edge diff=',ppoly_r_e(i0,1)-ppoly_r_e(i0-1,2)
435 problem_detected = .true.
438 if (problem_detected)
then 439 write(0,
'(a,1p9e24.16)')
'Polynomial coeffs:',ppoly_r_coefficients(i0,:)
440 write(0,
'(3(a,1pe24.16,x))')
'u_l=',u_l,
'u_c=',u_c,
'u_r=',u_r
441 write(0,
'(a4,10a24)')
'i0',
'h0(i0)',
'u0(i0)',
'left edge',
'right edge',
'Polynomial coefficients' 443 write(0,
'(i4,1p10e24.16)') n,h0(n),u0(n),ppoly_r_e(n,1),ppoly_r_e(n,2),ppoly_r_coefficients(n,:)
445 call mom_error(fatal,
'MOM_remapping, check_reconstructions_1d: '// &
446 'Edge values or polynomial coefficients were inconsistent!')
455 subroutine remap_via_sub_cells( n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, h1, method, &
456 force_bounds_in_subcell, u1, uh_err, ah_sub, aisub_src, aiss, aise )
457 integer,
intent(in) :: n0
458 real,
intent(in) :: h0(n0)
459 real,
intent(in) :: u0(n0)
460 real,
intent(in) :: ppoly0_E(n0,2)
461 real,
intent(in) :: ppoly0_coefficients(:,:)
462 integer,
intent(in) :: n1
463 real,
intent(in) :: h1(n1)
464 integer,
intent(in) :: method
465 logical,
intent(in) :: force_bounds_in_subcell
466 real,
intent(out) :: u1(n1)
467 real,
intent(out) :: uh_err
468 real,
optional,
intent(out) :: ah_sub(n0+n1+1)
469 integer,
optional,
intent(out) :: aisub_src(n0+n1+1)
470 integer,
optional,
intent(out) :: aiss(n0)
471 integer,
optional,
intent(out) :: aise(n0)
480 real,
dimension(n0+n1+1) :: h_sub
481 real,
dimension(n0+n1+1) :: uh_sub
482 real,
dimension(n0+n1+1) :: u_sub
483 integer,
dimension(n0+n1+1) :: isub_src
484 integer,
dimension(n0) :: isrc_start
485 integer,
dimension(n0) :: isrc_end
486 integer,
dimension(n0) :: isrc_max
487 real,
dimension(n0) :: h0_eff
488 real,
dimension(n0) :: u0_min
489 real,
dimension(n0) :: u0_max
490 integer,
dimension(n1) :: itgt_start
491 integer,
dimension(n1) :: itgt_end
493 real :: h0_supply, h1_supply
498 logical,
parameter :: force_bounds_in_target = .true.
499 logical,
parameter :: adjust_thickest_subcell = .true.
500 logical,
parameter :: debug_bounds = .false.
501 integer :: k, i0_last_thick_cell
502 real :: h0tot, h0err, h1tot, h1err, h2tot, h2err, u02_err
503 real :: u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, u2tot, u2err, u2min, u2max, u_orig
504 logical :: src_has_volume
505 logical :: tgt_has_volume
509 i0_last_thick_cell = 0
511 u0_min(i0) = min(ppoly0_e(i0,1), ppoly0_e(i0,2))
512 u0_max(i0) = max(ppoly0_e(i0,1), ppoly0_e(i0,2))
513 if (h0(i0)>0.) i0_last_thick_cell = i0
519 src_has_volume = .true.
520 tgt_has_volume = .true.
522 i_start0 = 1 ; i_start1 = 1
535 do i_sub = 2, n0+n1+1
539 dh = min(h0_supply, h1_supply)
544 dh0_eff = dh0_eff + min(dh, h0_supply)
553 if (dh >= dh_max)
then 560 if (h0_supply <= h1_supply .and. src_has_volume)
then 563 h1_supply = h1_supply - dh
566 isrc_start(i0) = i_start0
577 if (i0 < i0_last_thick_cell)
then 581 do while (h0_supply==0. .and. i0<i0_last_thick_cell)
596 src_has_volume = .false.
599 elseif (h0_supply >= h1_supply .and. tgt_has_volume)
then 602 h0_supply = h0_supply - dh
605 itgt_start(i1) = i_start1
617 tgt_has_volume = .false.
620 elseif (src_has_volume)
then 622 h_sub(i_sub) = h0_supply
624 isrc_start(i0) = i_start0
639 src_has_volume = .false.
641 elseif (tgt_has_volume)
then 643 h_sub(i_sub) = h1_supply
645 itgt_start(i1) = i_start1
654 tgt_has_volume = .false.
657 stop
'remap_via_sub_cells: THIS SHOULD NEVER HAPPEN!' 666 u_sub(1) = ppoly0_e(1,1)
679 dh0_eff = dh0_eff + dh
680 if (h0_eff(i0)>0.)
then 681 xb = dh0_eff / h0_eff(i0)
683 u_sub(i_sub) =
average_value_ppoly( n0, u0, ppoly0_e, ppoly0_coefficients, method, i0, xa, xb)
686 u_sub(i_sub) = u0(i0)
688 if (debug_bounds)
then 689 if (method<5 .and.(u_sub(i_sub)<u0_min(i0) .or. u_sub(i_sub)>u0_max(i0)))
then 690 write(0,*)
'Sub cell average is out of bounds',i_sub,
'method=',method
691 write(0,*)
'xa,xb: ',xa,xb
692 write(0,*)
'Edge values: ',ppoly0_e(i0,:),
'mean',u0(i0)
693 write(0,*)
'a_c: ',(u0(i0)-ppoly0_e(i0,1))+(u0(i0)-ppoly0_e(i0,2))
694 write(0,*)
'Polynomial coeffs: ',ppoly0_coefficients(i0,:)
695 write(0,*)
'Bounds min=',u0_min(i0),
'max=',u0_max(i0)
696 write(0,*)
'Average: ',u_sub(i_sub),
'rel to min=',u_sub(i_sub)-u0_min(i0),
'rel to max=',u_sub(i_sub)-u0_max(i0)
697 call mom_error( fatal,
'MOM_remapping, remap_via_sub_cells: '//&
698 'Sub-cell average is out of bounds!' )
701 if (force_bounds_in_subcell)
then 705 u_orig = u_sub(i_sub)
706 u_sub(i_sub) = max( u_sub(i_sub), u0_min(i0) )
707 u_sub(i_sub) = min( u_sub(i_sub), u0_max(i0) )
708 u02_err = u02_err + dh*abs( u_sub(i_sub) - u_orig )
710 uh_sub(i_sub) = dh * u_sub(i_sub)
712 if (isub_src(i_sub+1) /= i0)
then 721 u_sub(n0+n1+1) = ppoly0_e(n0,2)
722 uh_sub(n0+n1+1) = ppoly0_e(n0,2) * h_sub(n0+n1+1)
724 if (adjust_thickest_subcell)
then 728 do i0 = 1, i0_last_thick_cell
730 dh_max = h_sub(i_max)
731 if (dh_max > 0.)
then 734 do i_sub = isrc_start(i0), isrc_end(i0)
735 if (i_sub /= i_max) duh = duh + uh_sub(i_sub)
737 uh_sub(i_max) = u0(i0)*h0(i0) - duh
738 u02_err = u02_err + max( abs(uh_sub(i_max)), abs(u0(i0)*h0(i0)), abs(duh) )
746 if (h1(i1) > 0.)
then 748 i_sub = itgt_start(i1)
749 if (force_bounds_in_target)
then 753 do i_sub = itgt_start(i1), itgt_end(i1)
754 if (force_bounds_in_target)
then 755 u1min = min(u1min, u_sub(i_sub))
756 u1max = max(u1max, u_sub(i_sub))
758 dh = dh + h_sub(i_sub)
759 duh = duh + uh_sub(i_sub)
761 uh_err = uh_err + max(abs(duh),abs(uh_sub(i_sub)))*epsilon(duh)
765 uh_err = uh_err + abs(duh)*epsilon(duh)
766 if (force_bounds_in_target)
then 768 u1(i1) = max(u1min, min(u1max, u1(i1)))
770 uh_err = uh_err + dh*abs( u1(i1)-u_orig )
773 u1(i1) = u_sub(itgt_start(i1))
778 if (debug_bounds)
then 783 if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err+u02_err .and. abs(h1tot-h0tot)<h0err+h1err) &
784 .or. (abs(u2tot-u0tot)>u0err+u2err+u02_err .and. abs(h2tot-h0tot)<h0err+h2err) &
785 .or. (u1min<u0min .or. u1max>u0max) )
then 786 write(0,*)
'method = ',method
787 write(0,*)
'Source to sub-cells:' 788 write(0,*)
'H: h0tot=',h0tot,
'h2tot=',h2tot,
'dh=',h2tot-h0tot,
'h0err=',h0err,
'h2err=',h2err
789 if (abs(h2tot-h0tot)>h0err+h2err)
write(0,*)
'H non-conservation difference=',h2tot-h0tot,
'allowed err=',h0err+h2err,
' <-----!' 790 write(0,*)
'UH: u0tot=',u0tot,
'u2tot=',u2tot,
'duh=',u2tot-u0tot,
'u0err=',u0err,
'u2err=',u2err,
'adjustment err=',u02_err
791 if (abs(u2tot-u0tot)>u0err+u2err)
write(0,*)
'U non-conservation difference=',u2tot-u0tot,
'allowed err=',u0err+u2err,
' <-----!' 792 write(0,*)
'Sub-cells to target:' 793 write(0,*)
'H: h2tot=',h2tot,
'h1tot=',h1tot,
'dh=',h1tot-h2tot,
'h2err=',h2err,
'h1err=',h1err
794 if (abs(h1tot-h2tot)>h2err+h1err)
write(0,*)
'H non-conservation difference=',h1tot-h2tot,
'allowed err=',h2err+h1err,
' <-----!' 795 write(0,*)
'UH: u2tot=',u2tot,
'u1tot=',u1tot,
'duh=',u1tot-u2tot,
'u2err=',u2err,
'u1err=',u1err,
'uh_err=',uh_err
796 if (abs(u1tot-u2tot)>u2err+u1err)
write(0,*)
'U non-conservation difference=',u1tot-u2tot,
'allowed err=',u2err+u1err,
' <-----!' 797 write(0,*)
'Source to target:' 798 write(0,*)
'H: h0tot=',h0tot,
'h1tot=',h1tot,
'dh=',h1tot-h0tot,
'h0err=',h0err,
'h1err=',h1err
799 if (abs(h1tot-h0tot)>h0err+h1err)
write(0,*)
'H non-conservation difference=',h1tot-h0tot,
'allowed err=',h0err+h1err,
' <-----!' 800 write(0,*)
'UH: u0tot=',u0tot,
'u1tot=',u1tot,
'duh=',u1tot-u0tot,
'u0err=',u0err,
'u1err=',u1err,
'uh_err=',uh_err
801 if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err)
write(0,*)
'U non-conservation difference=',u1tot-u0tot,
'allowed err=',u0err+u1err+uh_err,
' <-----!' 802 write(0,*)
'U: u0min=',u0min,
'u1min=',u1min,
'u2min=',u2min
803 if (u1min<u0min)
write(0,*)
'U minimum overshoot=',u1min-u0min,
' <-----!' 804 if (u2min<u0min)
write(0,*)
'U2 minimum overshoot=',u2min-u0min,
' <-----!' 805 write(0,*)
'U: u0max=',u0max,
'u1max=',u1max,
'u2max=',u2max
806 if (u1max>u0max)
write(0,*)
'U maximum overshoot=',u1max-u0max,
' <-----!' 807 if (u2max>u0max)
write(0,*)
'U2 maximum overshoot=',u2max-u0max,
' <-----!' 808 write(0,
'(a3,6a24,2a3)')
'k',
'h0',
'left edge',
'u0',
'right edge',
'h1',
'u1',
'is',
'ie' 810 if (k<=min(n0,n1))
then 811 write(0,
'(i3,1p6e24.16,2i3)') k,h0(k),ppoly0_e(k,1),u0(k),ppoly0_e(k,2),h1(k),u1(k),itgt_start(k),itgt_end(k)
813 write(0,
'(i3,96x,1p2e24.16,2i3)') k,h1(k),u1(k),itgt_start(k),itgt_end(k)
815 write(0,
'(i3,1p4e24.16)') k,h0(k),ppoly0_e(k,1),u0(k),ppoly0_e(k,2)
818 write(0,
'(a3,2a24)')
'k',
'u0',
'Polynomial coefficients' 820 write(0,
'(i3,1p6e24.16)') k,u0(k),ppoly0_coefficients(k,:)
822 write(0,
'(a3,3a24,a3,2a24)')
'k',
'Sub-cell h',
'Sub-cell u',
'Sub-cell hu',
'i0',
'xa',
'xb' 828 dh0_eff = dh0_eff + dh
829 xb = dh0_eff / h0_eff(i0)
831 write(0,
'(i3,1p3e24.16,i3,1p2e24.16)') k,h_sub(k),u_sub(k),uh_sub(k),i0,xa,xb
833 if (isub_src(k+1) /= i0)
then 834 dh0_eff = 0.; xa = 0.
840 call mom_error( fatal,
'MOM_remapping, remap_via_sub_cells: '//&
841 'Remapping result is inconsistent!' )
847 uh_err = uh_err + u02_err
849 if (
present(ah_sub)) ah_sub(1:n0+n1+1) = h_sub(1:n0+n1+1)
850 if (
present(aisub_src)) aisub_src(1:n0+n1+1) = isub_src(1:n0+n1+1)
851 if (
present(aiss)) aiss(1:n0) = isrc_start(1:n0)
852 if (
present(aise)) aise(1:n0) = isrc_end(1:n0)
859 real function average_value_ppoly( n0, u0, ppoly0_E, ppoly0_coefficients, method, i0, xa, xb)
860 integer,
intent(in) :: n0
861 real,
intent(in) :: u0(:)
862 real,
intent(in) :: ppoly0_E(:,:)
863 real,
intent(in) :: ppoly0_coefficients(:,:)
864 integer,
intent(in) :: method
865 integer,
intent(in) :: i0
866 real,
intent(in) :: xa
867 real,
intent(in) :: xb
869 real :: u_ave, xa_2, xb_2, xa2pxb2, xapxb
870 real,
parameter :: r_3 = 1.0/3.0
872 real :: mx, a_L, a_R, u_c, Ya, Yb, my, xa2b2ab, Ya2b2ab, a_c
875 select case ( method )
880 ppoly0_coefficients(i0,1) &
881 + ppoly0_coefficients(i0,2) * 0.5 * ( xb + xa ) )
883 mx = 0.5 * ( xa + xb )
884 a_l = ppoly0_e(i0, 1)
885 a_r = ppoly0_e(i0, 2)
887 a_c = 0.5 * ( ( u_c - a_l ) + ( u_c - a_r ) )
890 xa2b2ab = (xa*xa+xb*xb)+xa*xb
891 u_ave = a_l + ( ( a_r - a_l ) * mx &
892 + a_c * ( 3. * ( xb + xa ) - 2.*xa2b2ab ) )
897 my = 0.5 * ( ya + yb )
898 ya2b2ab = (ya*ya+yb*yb)+ya*yb
899 u_ave = a_r + ( ( a_l - a_r ) * my &
900 + a_c * ( 3. * ( yb + ya ) - 2.*ya2b2ab ) )
905 xa2pxb2 = xa_2 + xb_2
908 ppoly0_coefficients(i0,1) &
909 + ( ppoly0_coefficients(i0,2) * 0.5 * ( xapxb ) &
910 + ( ppoly0_coefficients(i0,3) * r_3 * ( xa2pxb2 + xa*xb ) &
911 + ( ppoly0_coefficients(i0,4) * 0.25* ( xa2pxb2 * xapxb ) &
912 + ppoly0_coefficients(i0,5) * 0.2 * ( ( xb*xb_2 + xa*xa_2 ) * xapxb + xa_2*xb_2 ) ) ) ) )
914 call mom_error( fatal,
'The selected integration method is invalid' )
917 select case ( method )
919 u_ave = ppoly0_coefficients(i0,1)
923 a_l = ppoly0_e(i0, 1)
924 a_r = ppoly0_e(i0, 2)
927 u_ave = a_l + xa * ( a_r - a_l )
929 u_ave = a_r + ya * ( a_l - a_r )
935 a_l = ppoly0_e(i0, 1)
936 a_r = ppoly0_e(i0, 2)
938 a_c = 3. * ( ( u_c - a_l ) + ( u_c - a_r ) )
941 u_ave = a_l + xa * ( ( a_r - a_l ) + a_c * ya )
943 u_ave = a_r + ya * ( ( a_l - a_r ) + a_c * xa )
946 u_ave = ppoly0_coefficients(i0,1) &
947 + xa * ( ppoly0_coefficients(i0,2) &
948 + xa * ( ppoly0_coefficients(i0,3) &
949 + xa * ( ppoly0_coefficients(i0,4) &
950 + xa * ppoly0_coefficients(i0,5) ) ) )
952 call mom_error( fatal,
'The selected integration method is invalid' )
960 subroutine measure_input_bounds( n0, h0, u0, ppoly_E, h0tot, h0err, u0tot, u0err, u0min, u0max )
961 integer,
intent(in) :: n0
962 real,
dimension(n0),
intent(in) :: h0
963 real,
dimension(n0),
intent(in) :: u0
964 real,
dimension(n0,2),
intent(in) :: ppoly_E
965 real,
intent(out) :: h0tot
966 real,
intent(out) :: h0err
967 real,
intent(out) :: u0tot
968 real,
intent(out) :: u0err
969 real,
intent(out) :: u0min
970 real,
intent(out) :: u0max
978 u0tot = h0(1) * u0(1)
980 u0min = min( ppoly_e(1,1), ppoly_e(1,2) )
981 u0max = max( ppoly_e(1,1), ppoly_e(1,2) )
983 h0tot = h0tot + h0(k)
984 h0err = h0err + eps * max(h0tot, h0(k))
985 u0tot = u0tot + h0(k) * u0(k)
986 u0err = u0err + eps * max(abs(u0tot), abs(h0(k) * u0(k)))
987 u0min = min( u0min, ppoly_e(k,1), ppoly_e(k,2) )
988 u0max = max( u0max, ppoly_e(k,1), ppoly_e(k,2) )
995 integer,
intent(in) :: n1
996 real,
dimension(n1),
intent(in) :: h1
997 real,
dimension(n1),
intent(in) :: u1
998 real,
intent(out) :: h1tot
999 real,
intent(out) :: h1err
1000 real,
intent(out) :: u1tot
1001 real,
intent(out) :: u1err
1002 real,
intent(out) :: u1min
1003 real,
intent(out) :: u1max
1008 eps = epsilon(h1(1))
1011 u1tot = h1(1) * u1(1)
1016 h1tot = h1tot + h1(k)
1017 h1err = h1err + eps * max(h1tot, h1(k))
1018 u1tot = u1tot + h1(k) * u1(k)
1019 u1err = u1err + eps * max(abs(u1tot), abs(h1(k) * u1(k)))
1020 u1min = min(u1min, u1(k))
1021 u1max = max(u1max, u1(k))
1028 subroutine remapbyprojection( n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, h1, method, u1 )
1029 integer,
intent(in) :: n0
1030 real,
intent(in) :: h0(:)
1031 real,
intent(in) :: u0(:)
1032 real,
intent(in) :: ppoly0_E(:,:)
1033 real,
intent(in) :: ppoly0_coefficients(:,:)
1034 integer,
intent(in) :: n1
1035 real,
intent(in) :: h1(:)
1036 integer,
intent(in) :: method
1037 real,
intent(out) :: u1(:)
1053 xr = xl + h1(itarget)
1056 xl, xr, h1(itarget), u1(itarget), jstart, xstart )
1072 subroutine remapbydeltaz( n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, dx1, method, u1, h1 )
1073 integer,
intent(in) :: n0
1074 real,
intent(in) :: h0(:)
1075 real,
intent(in) :: u0(:)
1076 real,
intent(in) :: ppoly0_E(:,:)
1077 real,
intent(in) :: ppoly0_coefficients(:,:)
1078 integer,
intent(in) :: n1
1079 real,
intent(in) :: dx1(:)
1081 real,
intent(out) :: u1(:)
1082 real,
optional,
intent(out) :: h1(:)
1086 real :: xOld, hOld, uOld
1087 real :: xNew, hNew, h_err
1088 real :: uhNew, hFlux, uAve, fluxL, fluxR
1104 if (itarget == 0)
then 1108 elseif (itarget <= n0)
then 1109 xold = xold + h0(itarget)
1112 h_err = h_err + epsilon(hold) * max(hold, xold)
1117 xnew = xold + dx1(itarget+1)
1118 xl = min( xold, xnew )
1119 xr = max( xold, xnew )
1122 hflux = abs(dx1(itarget+1))
1124 xl, xr, hflux, uave, jstart, xstart )
1126 fluxr = dx1(itarget+1)*uave
1129 hnew = hold + ( dx1(itarget+1) - dx1(itarget) )
1130 hnew = max( 0., hnew )
1131 uhnew = ( uold * hold ) + ( fluxr - fluxl )
1133 u1(itarget) = uhnew / hnew
1137 if (
present(h1)) h1(itarget) = hnew
1147 xL, xR, hC, uAve, jStart, xStart )
1148 integer,
intent(in) :: n0
1149 real,
intent(in) :: h0(:)
1150 real,
intent(in) :: u0(:)
1151 real,
intent(in) :: ppoly0_E(:,:)
1152 real,
intent(in) :: ppoly0_coefficients(:,:)
1153 integer,
intent(in) :: method
1154 real,
intent(in) :: xL, xR
1155 real,
intent(in) :: hC
1156 real,
intent(out) :: uAve
1157 integer,
intent(inout) :: jStart
1159 real,
intent(inout) :: xStart
1168 real :: x0jLl, x0jLr
1169 real :: x0jRl, x0jRr
1172 real :: x0_2, x1_2, x02px12, x0px1
1173 real,
parameter :: r_3 = 1.0/3.0
1184 x0jlr = x0jll + h0(j)
1186 if ( ( xl >= x0jll ) .AND. ( xl <= x0jlr ) )
then 1203 'MOM_remapping, integrateReconOnInterval: '//&
1204 'The location of the left-most cell could not be found')
1216 if ( abs(xr - xl) == 0.0 )
then 1222 if ( h0(jl) == 0.0 )
then 1223 uave = 0.5 * ( ppoly0_e(jl,1) + ppoly0_e(jl,2) )
1228 select case ( method )
1230 uave = ppoly0_coefficients(jl,1)
1232 uave = ppoly0_coefficients(jl,1) &
1233 + xi0 * ppoly0_coefficients(jl,2)
1235 uave = ppoly0_coefficients(jl,1) &
1236 + xi0 * ( ppoly0_coefficients(jl,2) &
1237 + xi0 * ppoly0_coefficients(jl,3) )
1239 uave = ppoly0_coefficients(jl,1) &
1240 + xi0 * ( ppoly0_coefficients(jl,2) &
1241 + xi0 * ( ppoly0_coefficients(jl,3) &
1242 + xi0 * ( ppoly0_coefficients(jl,4) &
1243 + xi0 * ppoly0_coefficients(jl,5) ) ) )
1245 call mom_error( fatal,
'The selected integration method is invalid' )
1258 x0jrr = x0jrl + h0(j)
1260 if ( ( xr >= x0jrl ) .AND. ( xr <= x0jrr ) )
then 1269 if (xr>x0jrr) jr = n0
1275 if ( jl == jr )
then 1284 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__ 1285 xi0 = max( 0., min( 1., ( xl - x0jll ) / ( h0(jl) +
h_neglect ) ) )
1286 xi1 = max( 0., min( 1., ( xr - x0jll ) / ( h0(jl) +
h_neglect ) ) )
1288 xi0 = xl / h0(jl) - x0jll / ( h0(jl) +
h_neglect )
1289 xi1 = xr / h0(jl) - x0jll / ( h0(jl) +
h_neglect )
1292 hact = h0(jl) * ( xi1 - xi0 )
1297 select case ( method )
1299 q = ( xr - xl ) * ppoly0_coefficients(jl,1)
1301 q = ( xr - xl ) * ( &
1302 ppoly0_coefficients(jl,1) &
1303 + ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) )
1305 q = ( xr - xl ) * ( &
1306 ppoly0_coefficients(jl,1) &
1307 + ( ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) &
1308 + ppoly0_coefficients(jl,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1312 x02px12 = x0_2 + x1_2
1314 q = ( xr - xl ) * ( &
1315 ppoly0_coefficients(jl,1) &
1316 + ( ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) &
1317 + ( ppoly0_coefficients(jl,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1318 + ppoly0_coefficients(jl,4) * 0.25* ( x02px12 * x0px1 ) &
1319 + ppoly0_coefficients(jl,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1321 call mom_error( fatal,
'The selected integration method is invalid' )
1340 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__ 1341 xi0 = max( 0., min( 1., ( xl - x0jll ) / ( h0(jl) +
h_neglect ) ) )
1343 xi0 = (xl - x0jll) / ( h0(jl) +
h_neglect )
1347 hact = h0(jl) * ( xi1 - xi0 )
1349 select case ( method )
1351 q = q + ( x0jlr - xl ) * ppoly0_coefficients(jl,1)
1353 q = q + ( x0jlr - xl ) * ( &
1354 ppoly0_coefficients(jl,1) &
1355 + ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) )
1357 q = q + ( x0jlr - xl ) * ( &
1358 ppoly0_coefficients(jl,1) &
1359 + ( ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) &
1360 + ppoly0_coefficients(jl,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1364 x02px12 = x0_2 + x1_2
1366 q = q + ( x0jlr - xl ) * ( &
1367 ppoly0_coefficients(jl,1) &
1368 + ( ppoly0_coefficients(jl,2) * 0.5 * ( xi1 + xi0 ) &
1369 + ( ppoly0_coefficients(jl,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1370 + ppoly0_coefficients(jl,4) * 0.25* ( x02px12 * x0px1 ) &
1371 + ppoly0_coefficients(jl,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1373 call mom_error( fatal,
'The selected integration method is invalid' )
1377 if ( jr > (jl+1) )
then 1379 q = q + h0(k) * u0(k)
1386 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__ 1387 xi1 = max( 0., min( 1., ( xr - x0jrl ) / ( h0(jr) +
h_neglect ) ) )
1389 xi1 = (xr - x0jrl) / ( h0(jr) +
h_neglect )
1392 hact = hact + h0(jr) * ( xi1 - xi0 )
1394 select case ( method )
1396 q = q + ( xr - x0jrl ) * ppoly0_coefficients(jr,1)
1398 q = q + ( xr - x0jrl ) * ( &
1399 ppoly0_coefficients(jr,1) &
1400 + ppoly0_coefficients(jr,2) * 0.5 * ( xi1 + xi0 ) )
1402 q = q + ( xr - x0jrl ) * ( &
1403 ppoly0_coefficients(jr,1) &
1404 + ( ppoly0_coefficients(jr,2) * 0.5 * ( xi1 + xi0 ) &
1405 + ppoly0_coefficients(jr,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1409 x02px12 = x0_2 + x1_2
1411 q = q + ( xr - x0jrl ) * ( &
1412 ppoly0_coefficients(jr,1) &
1413 + ( ppoly0_coefficients(jr,2) * 0.5 * ( xi1 + xi0 ) &
1414 + ( ppoly0_coefficients(jr,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1415 + ppoly0_coefficients(jr,4) * 0.25* ( x02px12 * x0px1 ) &
1416 + ppoly0_coefficients(jr,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1418 call mom_error( fatal,
'The selected integration method is invalid' )
1424 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__ 1426 uave = ppoly0_coefficients(jl,1)
1440 integer,
intent(in) :: n1
1441 real,
dimension(:),
intent(in) :: h1
1442 integer,
intent(in) :: n2
1443 real,
dimension(:),
intent(in) :: h2
1444 real,
dimension(:),
intent(out) :: dx
1452 do k = 1, max(n1,n2)
1453 if (k <= n1) x1 = x1 + h1(k)
1464 check_reconstruction, check_remapping, force_bounds_in_subcell)
1467 character(len=*),
intent(in) :: remapping_scheme
1468 logical,
optional,
intent(in) :: boundary_extrapolation
1469 logical,
optional,
intent(in) :: check_reconstruction
1470 logical,
optional,
intent(in) :: check_remapping
1471 logical,
optional,
intent(in) :: force_bounds_in_subcell
1474 call remapping_set_param(cs, remapping_scheme=remapping_scheme, boundary_extrapolation=boundary_extrapolation, &
1475 check_reconstruction=check_reconstruction, check_remapping=check_remapping, &
1476 force_bounds_in_subcell=force_bounds_in_subcell)
1485 character(len=*),
intent(in) :: string
1510 call mom_error(fatal,
"setReconstructionType: "//&
1511 "Unrecognized choice for REMAPPING_SCHEME ("//trim(string)//
").")
1530 logical,
intent(in) :: verbose
1532 integer,
parameter :: n0 = 4, n1 = 3, n2 = 6
1533 real :: h0(n0), x0(n0+1), u0(n0)
1534 real :: h1(n1), x1(n1+1), u1(n1), hn1(n1), dx1(n1+1)
1535 real :: h2(n2), x2(n2+1), u2(n2), hn2(n2), dx2(n2+1)
1536 data u0 /9., 3., -3., -9./
1541 real,
allocatable,
dimension(:,:) :: ppoly0_E, ppoly0_S, ppoly0_coefficients
1544 logical :: thisTest, v
1548 write(*,*)
'==== MOM_remapping: remapping_unit_tests =================' 1554 err=x0(i)-0.75*
real(i-1)
1555 if (abs(err)>
real(i-1)*epsilon(err)) thistest = .true.
1557 if (thistest)
write(*,*)
'remapping_unit_tests: Failed buildGridFromH() 1' 1562 if (abs(err)>
real(i-1)*epsilon(err)) thistest = .true.
1564 if (thistest)
write(*,*)
'remapping_unit_tests: Failed buildGridFromH() 2' 1569 if (verbose)
write(*,*)
'h0 (test data)' 1570 if (verbose)
call dumpgrid(n0,h0,x0,u0)
1575 err=u1(i)-8.*(0.5*
real(1+n1)-
real(i))
1576 if (abs(err)>
real(n1-1)*epsilon(err)) thistest = .true.
1578 if (verbose)
write(*,*)
'h1 (by projection)' 1579 if (verbose)
call dumpgrid(n1,h1,x1,u1)
1580 if (thistest)
write(*,*)
'remapping_unit_tests: Failed remapping_core_w()' 1584 allocate(ppoly0_e(n0,2))
1585 allocate(ppoly0_s(n0,2))
1586 allocate(ppoly0_coefficients(n0,cs%degree+1))
1590 ppoly0_coefficients(:,:) = 0.0
1592 call edge_values_explicit_h4( n0, h0, u0, ppoly0_e )
1599 err=u1(i)-8.*(0.5*
real(1+n1)-
real(i))
1600 if (abs(err)>2.*epsilon(err)) thistest = .true.
1602 if (thistest)
write(*,*)
'remapping_unit_tests: Failed remapByProjection()' 1607 call remapbydeltaz( n0, h0, u0, ppoly0_e, ppoly0_coefficients, &
1608 n1, x1-x0(1:n1+1), &
1610 if (verbose)
write(*,*)
'h1 (by delta)' 1611 if (verbose)
call dumpgrid(n1,h1,x1,u1)
1614 err=u1(i)-8.*(0.5*
real(1+n1)-
real(i))
1615 if (abs(err)>2.*epsilon(err)) thistest = .true.
1617 if (thistest)
write(*,*)
'remapping_unit_tests: Failed remapByDeltaZ() 1' 1622 dx2(1:n0+1) = x2(1:n0+1) - x0
1623 dx2(n0+2:n2+1) = x2(n0+2:n2+1) - x0(n0+1)
1624 call remapbydeltaz( n0, h0, u0, ppoly0_e, ppoly0_coefficients, &
1627 if (verbose)
write(*,*)
'h2' 1628 if (verbose)
call dumpgrid(n2,h2,x2,u2)
1629 if (verbose)
write(*,*)
'hn2' 1630 if (verbose)
call dumpgrid(n2,hn2,x2,u2)
1633 err=u2(i)-8./2.*(0.5*
real(1+n2)-
real(i))
1634 if (abs(err)>2.*epsilon(err)) thistest = .true.
1636 if (thistest)
write(*,*)
'remapping_unit_tests: Failed remapByDeltaZ() 2' 1639 if (verbose)
write(*,*)
'Via sub-cells' 1643 if (verbose)
call dumpgrid(n2,h2,x2,u2)
1646 err=u2(i)-8./2.*(0.5*
real(1+n2)-
real(i))
1647 if (abs(err)>2.*epsilon(err)) thistest = .true.
1649 if (thistest)
write(*,*)
'remapping_unit_tests: Failed remap_via_sub_cells() 2' 1653 6, (/.125,.125,.125,.125,.125,.125/),
integration_ppm, .false., u2, err )
1654 if (verbose)
call dumpgrid(6,h2,x2,u2)
1658 if (verbose)
call dumpgrid(3,h2,x2,u2)
1662 write(*,*)
'===== MOM_remapping: new remapping_unit_tests ==================' 1664 deallocate(ppoly0_e, ppoly0_s, ppoly0_coefficients)
1665 allocate(ppoly0_coefficients(5,6))
1666 allocate(ppoly0_e(5,2))
1667 allocate(ppoly0_s(5,2))
1669 call pcm_reconstruction(3, (/1.,2.,4./), ppoly0_e(1:3,:), ppoly0_coefficients(1:3,:) )
1674 call plm_reconstruction(3, (/1.,1.,1./), (/1.,3.,5./), ppoly0_e(1:3,:), ppoly0_coefficients(1:3,:) )
1680 call plm_reconstruction(3, (/1.,1.,1./), (/1.,2.,7./), ppoly0_e(1:3,:), ppoly0_coefficients(1:3,:) )
1686 call plm_reconstruction(3, (/1.,1.,1./), (/1.,6.,7./), ppoly0_e(1:3,:), ppoly0_coefficients(1:3,:) )
1692 call plm_reconstruction(3, (/1.,2.,3./), (/1.,4.,9./), ppoly0_e(1:3,:), ppoly0_coefficients(1:3,:) )
1698 call edge_values_explicit_h4( 5, (/1.,1.,1.,1.,1./), (/1.,3.,5.,7.,9./), ppoly0_e )
1699 thistest =
test_answer(v, 5, ppoly0_e(:,1), (/0.,2.,4.,6.,8./),
'Line H4: left edges')
1700 thistest =
test_answer(v, 5, ppoly0_e(:,2), (/2.,4.,6.,8.,10./),
'Line H4: right edges')
1701 ppoly0_e(:,1) = (/0.,2.,4.,6.,8./)
1702 ppoly0_e(:,2) = (/2.,4.,6.,8.,10./)
1703 call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/1.,3.,5.,7.,9./), ppoly0_e(1:5,:), ppoly0_coefficients(1:5,:) )
1708 call edge_values_explicit_h4( 5, (/1.,1.,1.,1.,1./), (/1.,1.,7.,19.,37./), ppoly0_e )
1709 thistest =
test_answer(v, 5, ppoly0_e(:,1), (/3.,0.,3.,12.,27./),
'Parabola H4: left edges')
1710 thistest =
test_answer(v, 5, ppoly0_e(:,2), (/0.,3.,12.,27.,48./),
'Parabola H4: right edges')
1711 ppoly0_e(:,1) = (/0.,0.,3.,12.,27./)
1712 ppoly0_e(:,2) = (/0.,3.,12.,27.,48./)
1713 call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/0.,1.,7.,19.,37./), ppoly0_e(1:5,:), ppoly0_coefficients(1:5,:) )
1720 ppoly0_e(:,1) = (/0.,0.,6.,10.,15./)
1721 ppoly0_e(:,2) = (/0.,6.,12.,17.,15./)
1722 call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/0.,5.,7.,16.,15./), ppoly0_e(1:5,:), ppoly0_coefficients(1:5,:) )
1729 call plm_reconstruction(4, (/0.,1.,1.,0./), (/5.,4.,2.,1./), ppoly0_e(1:4,:), ppoly0_coefficients(1:4,:) )
1732 call remap_via_sub_cells( 4, (/0.,1.,1.,0./), (/5.,4.,2.,1./), ppoly0_e(1:4,:), ppoly0_coefficients(1:4,:), &
1736 deallocate(ppoly0_e, ppoly0_s, ppoly0_coefficients)
1743 logical function test_answer(verbose, n, u, u_true, label)
1744 logical,
intent(in) :: verbose
1745 integer,
intent(in) :: n
1746 real,
dimension(n),
intent(in) :: u
1747 real,
dimension(n),
intent(in) :: u_true
1748 character(len=*),
intent(in) :: label
1757 write(*,
'(a4,2a24,x,a)')
'k',
'Calculated value',
'Correct value',label
1759 if (u(k) /= u_true(k))
then 1760 write(*,
'(i4,1p2e24.16,a,1pe24.16,a)') k,u(k),u_true(k),
' err=',u(k)-u_true(k),
' < wrong' 1762 write(*,
'(i4,1p2e24.16)') k,u(k),u_true(k)
1771 integer,
intent(in) :: n
1772 real,
dimension(:),
intent(in) :: h
1773 real,
dimension(:),
intent(in) :: x
1774 real,
dimension(:),
intent(in) :: u
1776 write(*,
'("i=",20i10)') (i,i=1,n+1)
1777 write(*,
'("x=",20es10.2)') (x(i),i=1,n+1)
1778 write(*,
'("i=",5x,20i10)') (i,i=1,n)
1779 write(*,
'("h=",5x,20es10.2)') (h(i),i=1,n)
1780 write(*,
'("u=",5x,20es10.2)') (u(i),i=1,n)
subroutine, public edge_values_implicit_h6(N, h, u, edge_values)
subroutine, public edge_values_implicit_h4(N, h, u, edge_values)
character(len=3), public remappingdefaultscheme
Default remapping method.
subroutine buildgridfromh(nz, h, x)
Calculate edge coordinate x from cell width h.
subroutine measure_output_bounds(n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max)
Measure totals and bounds on destination grid.
subroutine remapbyprojection(n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, h1, method, u1)
Remaps column of values u0 on grid h0 to grid h1 by integrating over the projection of each h1 cell o...
integer, parameter remapping_ppm_h4
O(h^3) remapping scheme.
subroutine, public pqm_boundary_extrapolation_v1(N, h, u, ppoly_E, ppoly_S, ppoly_coefficients)
real function average_value_ppoly(n0, u0, ppoly0_E, ppoly0_coefficients, method, i0, xa, xb)
Returns the average value of a reconstruction within a single source cell, i0, between the non-dimens...
integer, parameter remapping_ppm_ih4
O(h^3) remapping scheme.
character(len=256), public remappingschemesdoc
Documentation for external callers.
Provides column-wise vertical remapping functions.
subroutine, public remapping_set_param(CS, remapping_scheme, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell)
Set parameters within remapping object.
subroutine setreconstructiontype(string, CS)
Changes the method of reconstruction Use this routine to parse a string parameter specifying the reco...
subroutine, public pqm_reconstruction(N, h, u, ppoly_E, ppoly_S, ppoly_coefficients)
subroutine, public pcm_reconstruction(N, u, ppoly_E, ppoly_coefficients)
Provides functions used with the Piecewise-Parabolic-Method in the vertical ALE algorithm.
subroutine check_reconstructions_1d(n0, h0, u0, deg, boundary_extrapolation, ppoly_r_coefficients, ppoly_r_E, ppoly_r_S)
Checks that edge values and reconstructions satisfy bounds.
logical function test_answer(verbose, n, u, u_true, label)
Returns true if any cell of u and u_true are not identical. Returns false otherwise.
logical, parameter old_algorithm
Use the old "broken" algorithm. This is a temporary measure to assist debugging until we delete the o...
logical function, public remapping_unit_tests(verbose)
Runs unit tests on remapping functions. Should only be called from a single/root thread Returns True ...
subroutine remap_via_sub_cells(n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, h1, method, force_bounds_in_subcell, u1, uh_err, ah_sub, aisub_src, aiss, aise)
Remaps column of n0 values u0 on grid h0 to grid h1 with n1 cells by calculating the n0+n1+1 sub-inte...
Container for remapping parameters.
subroutine, public remapping_core_w(CS, n0, h0, u0, n1, dx, u1)
Remaps column of values u0 on grid h0 to implied grid h1 where the interfaces of h1 differ from those...
integer, parameter remapping_plm
O(h^2) remapping scheme.
character(len=40) mdl
This module's name.
real, parameter h_neglect
A dimensional (H units) number that can be added to thicknesses in a denominator without changing the...
subroutine integraterecononinterval(n0, h0, u0, ppoly0_E, ppoly0_coefficients, method, xL, xR, hC, uAve, jStart, xStart)
Integrate the reconstructed column profile over a single cell.
subroutine, public ppm_reconstruction(N, h, u, ppoly_E, ppoly_coefficients)
Builds quadratic polynomials coefficients from cell mean and edge values.
subroutine build_reconstructions_1d(n0, h0, u0, remapping_scheme, deg, boundary_extrapolation, ppoly_r_coefficients, ppoly_r_E, ppoly_r_S, iMethod)
Creates polynomial reconstructions of u0 on the source grid h0.
character(len=len(input_string)) function, public uppercase(input_string)
subroutine, public edge_slopes_implicit_h3(N, h, u, edge_slopes)
subroutine, public dzfromh1h2(n1, h1, n2, h2, dx)
Calculates the change in interface positions based on h1 and h2.
subroutine, public edge_slopes_implicit_h5(N, h, u, edge_slopes)
subroutine, public plm_boundary_extrapolation(N, h, u, ppoly_E, ppoly_coefficients)
subroutine, public plm_reconstruction(N, h, u, ppoly_E, ppoly_coefficients)
integer, parameter remapping_pqm_ih6ih5
O(h^5) remapping scheme.
integer, parameter integration_pqm
Piecewise Quartic Method.
subroutine remapbydeltaz(n0, h0, u0, ppoly0_E, ppoly0_coefficients, n1, dx1, method, u1, h1)
Remaps column of values u0 on grid h0 to implied grid h1 where the interfaces of h1 differ from those...
subroutine, public end_remapping(CS)
Destrcutor for remapping control structure.
integer, parameter integration_ppm
Piecewise Parabolic Method.
subroutine, public remapping_core_h(CS, n0, h0, u0, n1, h1, u1)
Remaps column of values u0 on grid h0 to grid h1 assuming the top edge is aligned.
integer, parameter remapping_pqm_ih4ih3
O(h^4) remapping scheme.
integer, parameter remapping_pcm
O(h^1) remapping scheme.
integer, parameter integration_pcm
Piecewise Constant Method.
logical function ispossumerrsignificant(n1, sum1, n2, sum2)
Compare two summation estimates of positive data and judge if due to more than round-off. When two sums are calculated from different vectors that should add up to the same value, the results can differ by round off. The round off error can be bounded to be proportional to the number of operations. This function returns true if the difference between sum1 and sum2 is larger than than the estimated round off bound.
subroutine, public ppm_boundary_extrapolation(N, h, u, ppoly_E, ppoly_coefficients)
subroutine, public initialize_remapping(CS, remapping_scheme, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell)
Constructor for remapping control structure.
subroutine measure_input_bounds(n0, h0, u0, ppoly_E, h0tot, h0err, u0tot, u0err, u0min, u0max)
Measure totals and bounds on source grid.
subroutine, public mom_error(level, message, all_print)
integer, parameter integration_plm
Piecewise Linear Method.
subroutine, public edge_values_explicit_h4(N, h, u, edge_values)
subroutine dumpgrid(n, h, x, u)
Convenience function for printing grid to screen.