New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
rmap1d.h90 in vendors/PPR/src – NEMO

source: vendors/PPR/src/rmap1d.h90 @ 14212

Last change on this file since 14212 was 14212, checked in by jchanut, 4 years ago

#2222, update PPR library to circumvent issues found with cpp directives: remove timing instructions. Changed main module extension in .F90 to comply with FCM rules (no effect).

File size: 13.3 KB
Line 
1
2    !
3    ! This program may be freely redistributed under the
4    ! condition that the copyright notices (including this
5    ! entire header) are not removed, and no compensation
6    ! is received through use of the software.  Private,
7    ! research, and institutional use is free.  You may
8    ! distribute modified versions of this code UNDER THE
9    ! CONDITION THAT THIS CODE AND ANY MODIFICATIONS MADE
10    ! TO IT IN THE SAME FILE REMAIN UNDER COPYRIGHT OF THE
11    ! ORIGINAL AUTHOR, BOTH SOURCE AND OBJECT CODE ARE
12    ! MADE FREELY AVAILABLE WITHOUT CHARGE, AND CLEAR
13    ! NOTICE IS GIVEN OF THE MODIFICATIONS.  Distribution
14    ! of this code as part of a commercial system is
15    ! permissible ONLY BY DIRECT ARRANGEMENT WITH THE
16    ! AUTHOR.  (If you are not directly supplying this
17    ! code to a customer, and you are instead telling them
18    ! how they can obtain it for free, then you are not
19    ! required to make any arrangement with me.)
20    !
21    ! Disclaimer:  Neither I nor: Columbia University, the
22    ! National Aeronautics and Space Administration, nor
23    ! the Massachusetts Institute of Technology warrant
24    ! or certify this code in any way whatsoever.  This
25    ! code is provided "as-is" to be used at your own risk.
26    !
27    !
28
29    !   
30    ! RMAP1D.h90: high-order integral re-mapping operators.
31    !
32    ! Darren Engwirda
33    ! 31-Mar-2019
34    ! ​de2363 [at] columbia [dot] edu
35    !
36    !
37
38    subroutine rmap1d(npos,nnew,nvar,ndof,xpos, &
39        &             xnew,fdat,fnew,bclo,bcup, &
40        &             work,opts,tCPU)
41
42    !
43    ! NPOS  no. edges in old grid.
44    ! NNEW  no. edges in new grid.
45    ! NVAR  no. discrete variables to remap.
46    ! NDOF  no. degrees-of-freedom per cell.
47    ! XPOS  old grid edge positions. XPOS is a length NPOS
48    !       array .
49    ! XNEW  new grid edge positions. XNEW is a length NNEW
50    !       array .
51    ! FDAT  grid-cell moments on old grid. FNEW has SIZE =
52    !       NDOF-by-NVAR-by-NNEW-1 .
53    ! FNEW  grid-cell moments on new grid. FNEW has SIZE =
54    !       NDOF-by-NVAR-by-NNEW-1 .
55    ! BCLO  boundary condition at lower endpoint .
56    ! BCHI  boundary condition at upper endpoint .
57    ! WORK  method work-space. See RCON-WORK for details .
58    ! OPTS  method parameters. See RCON-OPTS for details .
59    ! TCPU  method tcpu-timer.
60    !
61
62        implicit none
63
64    !------------------------------------------- arguments !
65        integer, intent( in) :: npos,nnew
66        integer, intent( in) :: nvar,ndof
67        class(rmap_work), intent(inout):: work
68        class(rmap_opts), intent(inout):: opts
69        real*8 , intent( in) :: xpos(:)
70        real*8 , intent( in) :: xnew(:)
71        real*8 , intent( in) :: fdat(:,:,:)
72        real*8 , intent(out) :: fnew(:,:,:)
73        type (rcon_ends), intent(in) :: bclo(:)
74        type (rcon_ends), intent(in) :: bcup(:)
75        type (rmap_tics), &
76        &   intent(inout) , optional :: tCPU
77
78        real*8 , parameter :: RTOL = +1.d-14
79
80    !------------------------------------------- variables !
81        integer :: ipos
82        real*8  :: diff,spac,same,xtol
83        real*8  :: delx(1)
84        logical :: uniform
85       
86#       ifdef __PPR_TIMER__
87        integer(kind=8) :: ttic,ttoc,rate
88#       endif
89
90        if (ndof.lt.1) return
91        if (npos.lt.2) return
92        if (nnew.lt.2) return
93        if (nvar.lt.1) return
94
95    !------------- calc. grid-spacing and check uniformity !
96
97        same = (xpos(npos)&
98             -  xpos(  +1)) / (npos-1)
99
100        uniform = .true.
101             
102        xtol = same * RTOL
103
104        do  ipos = +1 , npos-1, +1
105       
106            spac = xpos(ipos+1) &
107        &        - xpos(ipos+0)
108 
109            diff = abs(spac - same)
110       
111            if (diff.gt.xtol) then
112           
113                uniform = .false.
114           
115            end if
116 
117            work% &
118        &   cell_spac(ipos) = spac
119       
120        end do
121
122       !uniform = .false.
123   
124    !------------- reconstruct FHAT over all cells in XPOS !
125
126        if (.not. uniform) then
127
128    !------------------------------------ variable spacing !
129        call rcon1d(npos,nvar,ndof, &
130        &           work%cell_spac, &
131        &           fdat,bclo,bcup, &
132        &           work%cell_func, &
133        &           work,opts,tCPU)
134
135        else
136       
137    !------------------------------------ constant spacing !
138        delx(1) = work%cell_spac(1)
139       
140        call rcon1d(npos,nvar,ndof, &
141        &           delx,    &
142        &           fdat,bclo,bcup, &
143        &           work%cell_func, &
144        &           work,opts,tCPU)
145       
146        end if
147
148    !------------- remap FDAT from XPOS to XNEW using FHAT !
149
150
151        select case(opts%cell_meth)
152        case(pcm_method)
153    !------------------------------------ 1st-order method !
154        call imap1d(npos,nnew,nvar, &
155        &           ndof,  +1,      &
156        &           xpos,xnew,      &
157        &           work%cell_func, &
158        &           fdat,fnew,xtol)
159   
160        case(plm_method)
161    !------------------------------------ 2nd-order method !
162        call imap1d(npos,nnew,nvar, &
163        &           ndof,  +2,      &
164        &           xpos,xnew,      &
165        &           work%cell_func, &
166        &           fdat,fnew,xtol)
167
168        case(ppm_method)
169    !------------------------------------ 3rd-order method !
170        call imap1d(npos,nnew,nvar, &
171        &           ndof,  +3,      &
172        &           xpos,xnew,      &
173        &           work%cell_func, &
174        &           fdat,fnew,xtol)
175
176        case(pqm_method)
177    !------------------------------------ 5th-order method !
178        call imap1d(npos,nnew,nvar, &
179        &           ndof,  +5,      &
180        &           xpos,xnew,      &
181        &           work%cell_func, &
182        &           fdat,fnew,xtol)
183       
184        end select
185       
186       
187        return
188   
189    end  subroutine
190   
191    !------------ IMAP1D: 1-dimensional degree-k remapping !
192   
193    pure subroutine imap1d(npos,nnew,nvar,ndof, &
194        &                  mdof,xpos,xnew,fhat, &
195        &                  fdat,fnew,XTOL)
196
197    !
198    ! NPOS  no. edges in old grid.
199    ! NNEW  no. edges in new grid.
200    ! NVAR  no. discrete variables to remap.
201    ! NDOF  no. degrees-of-freedom per cell.
202    ! MDOF  no. degrees-of-freedom per FHAT.   
203    ! XPOS  old grid edge positions. XPOS is a length NPOS
204    !       array .
205    ! XNEW  new grid edge positions. XNEW is a length NNEW
206    !       array .
207    ! FHAT  reconstruction over old grid. FHAT has SIZE =
208    !       MDOF-by-NVAR-by-NPOS-1 .
209    ! FDAT  grid-cell moments on old grid. FDAT has SIZE =
210    !       NDOF-by-NVAR-by-NPOS-1 .
211    ! FNEW  grid-cell moments on new grid. FNEW has SIZE =
212    !       NDOF-by-NVAR-by-NNEW-1 .
213    ! XTOL  min. grid-cell thickness .
214    !
215
216        implicit none   
217   
218    !------------------------------------------- arguments !
219        integer, intent( in) :: npos,nnew
220        integer, intent( in) :: nvar
221        integer, intent( in) :: ndof,mdof
222        real*8 , intent( in) :: xpos(:)
223        real*8 , intent( in) :: xnew(:)
224        real*8 , intent( in) :: fhat(:,:,:)
225        real*8 , intent( in) :: fdat(:,:,:)       
226        real*8 , intent(out) :: fnew(:,:,:)
227        real*8 , intent( in) :: XTOL
228       
229    !------------------------------------------- variables !
230        integer :: kpos,ipos,ivar,idof
231        integer :: pos0,pos1,vmin,vmax
232        real*8  :: xmid,xhat,khat,stmp
233        real*8  :: xxlo,xxhi,sslo,sshi,intf
234        real*8  :: vvlo(  +1:+5)
235        real*8  :: vvhi(  +1:+5)       
236        real*8  :: ivec(  +1:+5)
237        real*8  :: sdat(  +1:nvar)
238        real*8  :: snew(  +1:nvar)
239        real*8  :: serr(  +1:nvar)
240        integer :: kmin(  +1:nvar)
241        integer :: kmax(  +1:nvar)
242       
243        integer, parameter :: INTB = -1   ! integral basis 
244
245    !------------------------------------- initializations !
246
247        vvlo(+1:+5) = 0.d0
248        vvhi(+1:+5) = 0.d0
249
250    !------------- remap FDAT from XPOS to XNEW using FHAT !
251
252        kmin = +1 ; kmax = +1
253        pos0 = +1 ; pos1 = +1
254
255        do  kpos = +1, nnew-1
256
257    !------ first cell in XPOS overlapping with XNEW(KPOS) !
258
259            pos1 = max(pos1,1)
260
261            do  pos0 = pos1, npos-1
262
263                if (xpos(pos0+1)&
264            &       .gt. xnew(kpos+0)) exit
265
266            end do
267
268    !------ final cell in XPOS overlapping with XNEW(KPOS) !
269   
270            do  pos1 = pos0, npos-1
271
272                if (xpos(pos1+0)&
273            &       .ge. xnew(kpos+1)) exit
274   
275            end do
276           
277            pos1 = pos1 - 1
278
279    !------------- integrate FHAT across overlapping cells !
280
281            khat = xnew(kpos+1) &
282            &    - xnew(kpos+0)
283            khat = max (khat , XTOL)
284
285            do  idof = +1,ndof
286            do  ivar = +1,nvar
287               
288                fnew(idof,ivar,kpos) = 0.d0
289           
290            end do
291            end do
292
293            do  ipos = pos0, pos1
294
295    !------------------------------- integration endpoints !
296   
297                xxlo = max (xpos(ipos+0) , &
298            &               xnew(kpos+0))
299                xxhi = min (xpos(ipos+1) , &
300            &               xnew(kpos+1))
301
302    !------------------------------- local endpoint coords !
303   
304                xmid = xpos(ipos+1) * .5d0 &
305            &        + xpos(ipos+0) * .5d0   
306                xhat = xpos(ipos+1) * .5d0 &
307            &        - xpos(ipos+0) * .5d0
308     
309                sslo = &
310            &  (xxlo-xmid) / max(xhat,XTOL)
311                sshi = &
312            &  (xxhi-xmid) / max(xhat,XTOL)
313
314    !------------------------------- integral basis vector !
315   
316                call bfun1d(INTB,mdof, &
317                            sslo,vvlo)
318                call bfun1d(INTB,mdof, &
319                            sshi,vvhi)
320               
321                ivec =  vvhi - vvlo
322
323    !--------- integrate FHAT across the overlap XXLO:XXHI !
324   
325                do  ivar = +1, nvar
326   
327                intf =  dot_product (  &
328            &       ivec(+1:mdof),  &
329            &   fhat(+1:mdof,ivar,ipos-0) )
330
331                intf =  intf * xhat
332       
333    !--------- accumulate integral contributions from IPOS !
334   
335                fnew(  +1,ivar,kpos) = &
336            &   fnew(  +1,ivar,kpos) + intf
337
338                end do
339
340            end do
341
342    !------------------------------- finalise KPOS profile !
343
344            do  ivar = +1, nvar
345
346                fnew(  +1,ivar,kpos) = &
347            &   fnew(  +1,ivar,kpos) / khat
348
349    !--------- keep track of MIN/MAX for defect correction !
350
351                vmax =    kmax(ivar)
352                vmin =    kmin(ivar)
353
354                if(fnew(1,ivar,kpos) &
355            &  .gt.fnew(1,ivar,vmax) ) then
356               
357                kmax(ivar) =   kpos
358           
359                else &
360            &   if(fnew(1,ivar,kpos) &
361            &  .lt.fnew(1,ivar,vmin) ) then
362               
363                kmin(ivar) =   kpos
364           
365                end if
366
367            end do
368
369        end do
370
371    !--------- defect corrections: Kahan/Babuska/Neumaier. !
372
373    !   Carefully compute column sums, leading to a defect
374    !   wrt. column-wise conservation. Use KBN approach to
375    !   account for FP roundoff.
376
377        sdat = 0.d0; serr = 0.d0
378        do  ipos = +1, npos-1
379        do  ivar = +1, nvar-0
380       
381    !------------------------------- integrate old profile !
382
383            xhat = xpos(ipos+1) &
384        &        - xpos(ipos+0)
385
386            intf = xhat*fdat(1,ivar,ipos)
387           
388            stmp = sdat(ivar) + intf
389
390            if (abs(sdat(ivar)) &
391        &           .ge. abs(intf)) then
392
393            serr(ivar) = &
394        &   serr(ivar) + ((sdat(ivar)-stmp)+intf)
395
396            else
397
398            serr(ivar) = &
399        &   serr(ivar) + ((intf-stmp)+sdat(ivar))
400
401            end if
402
403            sdat(ivar) = stmp
404       
405        end do
406        end do
407
408        sdat =  sdat + serr
409
410        snew = 0.d0; serr = 0.d0
411        do  ipos = +1, nnew-1
412        do  ivar = +1, nvar-0
413       
414    !------------------------------- integrate new profile !
415
416            khat = xnew(ipos+1) &
417        &        - xnew(ipos+0)
418
419            intf = khat*fnew(1,ivar,ipos)
420           
421            stmp = snew(ivar) + intf
422
423            if (abs(snew(ivar)) &
424        &           .ge. abs(intf)) then
425
426            serr(ivar) = &
427        &   serr(ivar) + ((snew(ivar)-stmp)+intf)
428
429            else
430
431            serr(ivar) = &
432        &   serr(ivar) + ((intf-stmp)+snew(ivar))
433
434            end if
435
436            snew(ivar) = stmp
437       
438        end do
439        end do
440
441        snew =  snew + serr
442        serr =  sdat - snew
443
444    !--------- defect corrections: nudge away from extrema !
445
446    !   Add a correction to remapped state to impose exact
447    !   conservation. Via sign(correction), nudge min/max.
448    !   cell means, such that monotonicity is not violated
449    !   near extrema...
450   
451        do  ivar = +1, nvar-0
452
453            if (serr(ivar) .gt. 0.d0) then
454
455                vmin = kmin(ivar)
456
457                fnew(1,ivar,vmin) = &
458        &       fnew(1,ivar,vmin) + &
459        &  serr(ivar)/(xnew(vmin+1)-xnew(vmin+0))
460
461            else &
462        &   if (serr(ivar) .lt. 0.d0) then
463
464                vmax = kmax(ivar)
465
466                fnew(1,ivar,vmax) = &
467        &       fnew(1,ivar,vmax) + &
468        &  serr(ivar)/(xnew(vmax+1)-xnew(vmax+0))
469
470            end if
471
472        end do
473       
474    !------------------------------- new profile now final !
475
476        return
477   
478    end  subroutine
479
480
481
Note: See TracBrowser for help on using the repository browser.