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.
ppr_1d.F90 in vendors/PPR/src – NEMO

source: vendors/PPR/src/ppr_1d.F90 @ 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: 9.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    module ppr_1d
30
31    !   
32    ! PPR-1D.F90: 1-d piecewise polynomial reconstructions.
33    !
34    ! Darren Engwirda
35    ! 31-Mar-2019
36    ! de2363 [at] columbia [dot] edu
37    !
38    !
39
40    implicit none
41
42    !------------------------------------ method selection !
43       
44    integer, parameter :: p1e_method = +100
45    integer, parameter :: p3e_method = +101
46    integer, parameter :: p5e_method = +102
47
48    integer, parameter :: pcm_method = +200
49    integer, parameter :: plm_method = +201
50    integer, parameter :: ppm_method = +202
51    integer, parameter :: pqm_method = +203
52
53    integer, parameter :: null_limit = +300
54    integer, parameter :: mono_limit = +301
55    integer, parameter :: weno_limit = +302
56
57    integer, parameter :: bcon_loose = +400
58    integer, parameter :: bcon_value = +401
59    integer, parameter :: bcon_slope = +402
60 
61    type rmap_tics
62    !------------------------------- tCPU timer for RCON1D !
63        integer :: rmap_time
64        integer :: edge_time
65        integer :: cell_time
66        integer :: oscl_time
67    end type rmap_tics
68
69    type rcon_opts
70    !------------------------------- parameters for RCON1D !
71        integer :: edge_meth
72        integer :: cell_meth
73        integer :: cell_lims
74        integer :: wall_lims
75    end type rcon_opts
76
77    type rcon_ends
78    !------------------------------- end-conditions struct !
79        integer :: bcopt
80        real*8  :: value
81        real*8  :: slope
82    end type rcon_ends
83
84    type rcon_work
85    !------------------------------- work-space for RCON1D !
86        real*8, allocatable  :: edge_func(:,:)
87        real*8, allocatable  :: edge_dfdx(:,:)
88        real*8, allocatable  :: cell_oscl(:,:,:)
89    contains
90        procedure :: init => init_rcon_work
91        procedure :: free => free_rcon_work
92    end type rcon_work
93
94    type, extends(rcon_opts) :: rmap_opts
95    !------------------------------- parameters for RMAP1D !
96    end type rmap_opts
97
98    type, extends(rcon_work) :: rmap_work
99    !------------------------------- work-space for RMAP1D !
100        real*8, allocatable  :: cell_spac(:)
101        real*8, allocatable  :: cell_func(:,:,:)
102    contains
103        procedure :: init => init_rmap_work
104        procedure :: free => free_rmap_work
105    end type rmap_work
106
107    contains
108 
109    !------------------------------------------------------!
110    ! INIT-RCON-WORK: init. work-space for RCON1D.         !
111    !------------------------------------------------------!
112   
113    subroutine init_rcon_work(this,npos,nvar,opts)
114
115    !
116    ! THIS  work-space structure for RCON1D .
117    ! NPOS  no. edges over grid.
118    ! NVAR  no. state variables.
119    ! OPTS  parameters structure for RCON1D .
120    !
121
122        implicit none
123
124    !------------------------------------------- arguments !
125        class(rcon_work) , intent(inout) :: this
126        integer, intent(in):: npos
127        integer, intent(in):: nvar
128        class(rcon_opts) , optional      :: opts
129
130    !------------------------------------------- variables !
131        integer :: okay
132
133        allocate(this% &
134        &   edge_func(  nvar,npos), &
135        &        this% &
136        &   edge_dfdx(  nvar,npos), &
137        &        this% &
138        &   cell_oscl(2,nvar,npos), &
139        &   stat=okay)
140       
141    end subroutine
142
143    !------------------------------------------------------!
144    ! INIT-RMAP-WORK: init. work-space for RMAP1D.         !
145    !------------------------------------------------------!
146   
147    subroutine init_rmap_work(this,npos,nvar,opts)
148
149    !
150    ! THIS  work-space structure for RMAP1D .
151    ! NPOS  no. edges over grid.
152    ! NVAR  no. state variables.
153    ! OPTS  parameters structure for RMAP1D .
154    !
155
156        implicit none
157
158    !------------------------------------------- arguments !
159        class(rmap_work) , intent(inout) :: this
160        integer, intent(in) :: npos
161        integer, intent(in) :: nvar
162        class(rcon_opts) , optional      :: opts
163
164    !------------------------------------------- variables !
165        integer :: okay,ndof
166
167        ndof = ndof1d(opts%cell_meth)
168
169        allocate(this% &
170        &   edge_func(  nvar,npos), &
171        &        this% &
172        &   edge_dfdx(  nvar,npos), &
173        &        this% &
174        &   cell_oscl(2,nvar,npos), &
175        &        this% &
176        &   cell_spac(       npos), &
177        &        this% &
178        &   cell_func(ndof,nvar,npos) , &
179        &   stat=okay)
180 
181    end subroutine
182
183    !------------------------------------------------------!
184    ! FREE-RCON-WORK: free work-space for RCON1D .         !
185    !------------------------------------------------------!
186   
187    subroutine free_rcon_work(this)
188
189        implicit none
190
191    !------------------------------------------- arguments !
192        class(rcon_work), intent(inout) :: this
193
194        deallocate(this%edge_func, &
195        &          this%edge_dfdx, &
196        &          this%cell_oscl)
197     
198    end subroutine
199
200    !------------------------------------------------------!
201    ! FREE-RMAP-WORK: free work-space for RMAP1D .         !
202    !------------------------------------------------------!
203   
204    subroutine free_rmap_work(this)
205
206        implicit none
207
208    !------------------------------------------- arguments !
209        class(rmap_work), intent(inout) :: this
210
211
212        deallocate(this%edge_func, &
213        &          this%edge_dfdx, &
214        &          this%cell_oscl, &
215        &          this%cell_func, &
216        &          this%cell_spac)
217
218    end subroutine
219
220    !------------------------------------------------------!
221    ! NDOF1D : no. degrees-of-freedom per polynomial .     !
222    !------------------------------------------------------!
223   
224    pure function ndof1d(meth) result(rdof)
225
226        implicit none
227
228    !------------------------------------------- arguments !
229        integer, intent( in) :: meth
230
231    !------------------------------------------- variables !
232        integer  :: rdof
233
234        select case(meth)
235    !-------------------------------- edge reconstructions !
236        case (p1e_method)
237            rdof = +2
238        case (p3e_method)
239            rdof = +4
240        case (p5e_method)
241            rdof = +6
242    !-------------------------------- cell reconstructions !
243        case (pcm_method)
244            rdof = +1
245        case (plm_method)
246            rdof = +2
247        case (ppm_method)
248            rdof = +3
249        case (pqm_method)
250            rdof = +5
251       
252        case default 
253            rdof = +0
254
255        end select
256       
257    end function  ndof1d
258
259    !------------------------------------------------------!
260    ! BFUN1D : one-dimensional poly. basis-functions .     !
261    !------------------------------------------------------!
262       
263#       include "bfun1d.h90"
264
265    !------------------------------------------------------!
266    ! UTIL1D : one-dimensional grid manip. utilities .     !
267    !------------------------------------------------------!
268   
269#       include "util1d.h90"
270
271    !------------------------------------------------------!
272    ! WENO1D : "essentially" non-oscillatory limiter .     !
273    !------------------------------------------------------!
274
275#       include "weno1d.h90"
276
277#       include "oscl1d.h90"
278
279    !------------------------------------------------------!
280    ! RCON1D : one-dimensional poly. reconstructions .     !
281    !------------------------------------------------------!
282
283#       include "rcon1d.h90"
284
285#       include "inv.h90"
286
287#       include "pbc.h90"
288#       include "p1e.h90"
289#       include "p3e.h90"
290#       include "p5e.h90"
291
292#       include "root1d.h90"
293   
294#       include "pcm.h90"
295#       include "plm.h90"
296#       include "ppm.h90"
297#       include "pqm.h90"
298   
299    !------------------------------------------------------!
300    ! RMAP1D : one-dimensional conservative "re-map" .     !
301    !------------------------------------------------------!
302       
303#       include "rmap1d.h90"
304
305    !------------------------------------------------------!
306    ! FFSL1D : one-dimensional FFSL scalar transport .     !
307    !------------------------------------------------------!
308   
309#       include "ffsl1d.h90"
310
311
312    !------------------------------------------ end ppr_1d !
313     
314       
315    end module
316   
317   
318   
Note: See TracBrowser for help on using the repository browser.