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.
rcon1d.h90 in NEMO/branches/2020/dev_r13312_AGRIF-03-04_jchanut_vinterp_tstep/ext/PPR/src – NEMO

source: NEMO/branches/2020/dev_r13312_AGRIF-03-04_jchanut_vinterp_tstep/ext/PPR/src/rcon1d.h90 @ 13926

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

#2222, add Piecewise Polynomial Reconstruction library

File size: 6.4 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    ! RCON1D.f90: conservative, polynomial reconstructions.
31    !
32    ! Darren Engwirda
33    ! 07-Sep-2016
34    ! de2363 [at] columbia [dot] edu
35    !
36    !
37   
38    subroutine rcon1d(npos,nvar,ndof,delx,fdat, &
39        &             bclo,bchi,fhat,work,opts, &
40        &             tCPU)
41
42    !
43    ! NPOS  no. edges over grid.
44    ! NVAR  no. state variables.
45    ! NDOF  no. degrees-of-freedom per grid-cell.
46    ! DELX  grid-cell spacing array. LENGTH(DELX) == +1 if
47    !       spacing is uniform .
48    ! FDAT  grid-cell moments array. FDAT is an array with
49    !       SIZE = NDOF-by-NVAR-by-NPOS-1 .
50    ! BCLO  boundary condition at lower endpoint.
51    ! BCHI  boundary condition at upper endpoint.
52    ! FHAT  grid-cell re-con. array. FHAT is an array with
53    !       SIZE = MDOF-by-NVAR-by-NPOS-1 .
54    ! WORK  method work-space. See RCON-WORK for details .
55    ! OPTS  method parameters. See RCON-OPTS for details .
56    ! TCPU  method tcpu-timer.
57    !
58
59        implicit none
60
61    !------------------------------------------- arguments !
62        integer, intent( in) :: npos,nvar,ndof
63        class(rcon_work), intent(inout):: work
64        class(rcon_opts), intent(in)   :: opts
65        real*8 , intent( in) :: delx(:)
66        real*8 , intent(out) :: fhat(:,:,:)
67        real*8 , intent( in) :: fdat(:,:,:)
68        type (rcon_ends), intent(in) :: bclo(:)
69        type (rcon_ends), intent(in) :: bchi(:)
70        type (rmap_tics), &
71        &   intent(inout) , optional :: tCPU
72
73    !------------------------------------------- variables !
74        integer :: halo,ipos
75        real*8  :: dmin,dmid
76
77#       ifdef __PPR_TIMER__
78        integer(kind=8) :: ttic,ttoc,rate
79#       endif
80
81        if (ndof.lt.1) return
82        if (npos.lt.2) return
83        if (nvar.lt.1) return
84     
85    !-------------------------- compute min grid-tolerance !
86     
87        dmid = delx(1)
88     
89        if (size(delx).gt.+1) then
90       
91            do  ipos = 2, npos-1
92                dmid = &
93            &   dmid + delx (ipos)
94            end do
95       
96            dmid = dmid /(npos-1)
97       
98        end if
99
100        dmin = +1.0d-14 * dmid
101       
102    !-------------------------- compute edge values/slopes !
103
104        __TIC__
105
106        if ( (opts%cell_meth.eq.ppm_method) &
107    &  .or.  (opts%cell_meth.eq.pqm_method) ) then
108
109        select case (opts%edge_meth)
110            case(p1e_method)
111    !------------------------------------ 2nd-order method !
112            halo = +1
113            call p1e(npos,nvar,ndof, &
114            &        delx,fdat,      &
115            &        bclo,bchi,      &
116            &        work%edge_func, &
117            &        work%edge_dfdx, &
118            &        opts,dmin)
119
120            case(p3e_method)
121    !------------------------------------ 4th-order method !           
122            halo = +2
123            call p3e(npos,nvar,ndof, &
124            &        delx,fdat,      &
125            &        bclo,bchi,      &
126            &        work%edge_func, &
127            &        work%edge_dfdx, &
128            &        opts,dmin)
129
130            case(p5e_method)
131    !------------------------------------ 6th-order method !           
132            halo = +3
133            call p5e(npos,nvar,ndof, &
134            &        delx,fdat,      &
135            &        bclo,bchi,      &
136            &        work%edge_func, &
137            &        work%edge_dfdx, &
138            &        opts,dmin)
139
140        end select
141
142        end if
143       
144        __TOC__(tCPU,edge_time)
145
146    !-------------------------- compute oscil. derivatives !
147
148        __TIC__
149
150        if (opts%cell_lims.eq.weno_limit) then
151   
152            call oscli(npos,nvar,ndof, &
153            &          delx,fdat, &
154            &          work%cell_oscl, &
155            &          dmin)
156     
157        end if
158   
159        __TOC__(time,oscl_time)
160
161    !-------------------------- compute grid-cell profiles !
162
163        __TIC__
164
165        select case (opts%cell_meth)
166            case(pcm_method)
167    !------------------------------------ 1st-order method !
168            call pcm(npos,nvar,ndof, &
169            &        fdat,fhat)
170
171            case(plm_method)
172    !------------------------------------ 2nd-order method !
173            call plm(npos,nvar,ndof, &
174            &        delx,fdat,fhat, &
175            &        dmin,&
176            &        opts%cell_lims)
177
178            case(ppm_method)
179    !------------------------------------ 3rd-order method !
180            call ppm(npos,nvar,ndof, &
181            &        delx,fdat,fhat, &
182            &        work%edge_func, &
183            &        work%cell_oscl, &
184            &        dmin,&
185            &        opts%cell_lims, &
186            &        opts%wall_lims, &
187            &        halo )
188
189            case(pqm_method)
190    !------------------------------------ 5th-order method !
191            call pqm(npos,nvar,ndof, &
192            &        delx,fdat,fhat, &
193            &        work%edge_func, &
194            &        work%edge_dfdx, &
195            &        work%cell_oscl, &
196            &        dmin,&
197            &        opts%cell_lims, &
198            &        opts%wall_lims, &
199            &        halo )
200
201        end select
202
203        __TOC__(tCPU,cell_time)
204
205    end subroutine
206   
207   
208   
Note: See TracBrowser for help on using the repository browser.