Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
ab_time_scheme.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2!
3! The UChicago Argonne, LLC as Operator of Argonne National
4! Laboratory holds copyright in the Software. The copyright holder
5! reserves all rights except those expressly granted to licensees,
6! and U.S. Government license rights.
7!
8! Redistribution and use in source and binary forms, with or without
9! modification, are permitted provided that the following conditions
10! are met:
11!
12! 1. Redistributions of source code must retain the above copyright
13! notice, this list of conditions and the disclaimer below.
14!
15! 2. Redistributions in binary form must reproduce the above copyright
16! notice, this list of conditions and the disclaimer (as noted below)
17! in the documentation and/or other materials provided with the
18! distribution.
19!
20! 3. Neither the name of ANL nor the names of its contributors
21! may be used to endorse or promote products derived from this software
22! without specific prior written permission.
23!
24! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36!
37! Additional BSD Notice
38! ---------------------
39! 1. This notice is required to be provided under our contract with
40! the U.S. Department of Energy (DOE). This work was produced at
41! Argonne National Laboratory under Contract
42! No. DE-AC02-06CH11357 with the DOE.
43!
44! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45! LLC nor any of their employees, makes any warranty,
46! express or implied, or assumes any liability or responsibility for the
47! accuracy, completeness, or usefulness of any information, apparatus,
48! product, or process disclosed, or represents that its use would not
49! infringe privately-owned rights.
50!
51! 3. Also, reference herein to any specific commercial products, process,
52! or services by trade name, trademark, manufacturer or otherwise does
53! not necessarily constitute or imply its endorsement, recommendation,
54! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55! The views and opinions of authors expressed
56! herein do not necessarily state or reflect those of the United States
57! Government or UCHICAGO ARGONNE, LLC, and shall
58! not be used for advertising or product endorsement purposes.
59!
62 use neko_config
63 use num_types, only : rp
64 use time_scheme, only: time_scheme_t
65 use math, only : rzero
66 use utils, only : neko_error
67 implicit none
68 private
69
76 type, public, extends(time_scheme_t) :: ab_time_scheme_t
77 contains
80 end type ab_time_scheme_t
81
82contains
83
87 subroutine ab_time_scheme_compute_coeffs(coeffs, dt, order)
88 implicit none
89 real(kind=rp), intent(out) :: coeffs(4)
90 real(kind=rp), intent(in) :: dt(10)
91 integer, intent(in) :: order
92 real(kind=rp) dta, dtb, dtc, dtd, dte, dts
93
94 call rzero(coeffs, 4)
95
96 select case (order)
97 case (1)
98 coeffs(1) = 1.0_rp
99 case (2)
100 coeffs(2) = -0.5_rp * dt(1) / dt(2)
101 coeffs(1) = 1.0_rp - coeffs(2)
102 case (3)
103 dts = dt(2) + dt(3)
104 dta = dt(1) / dt(2)
105 dtb = dt(2) / dt(3)
106 dtc = dt(1) / dt(3)
107 dtd = dts / dt(2)
108 dte = dt(1) / dts
109 coeffs(3) = dte*( 0.5_rp*dtb + dtc/3.0_rp )
110 coeffs(2) = -0.5_rp * dta - coeffs(3) * dtd
111 coeffs(1) = 1.0_rp - coeffs(2) - coeffs(3)
112 case default
113 call neko_error("The order of the AB time scheme must be 1 to 3.")
114 end select
115
116 end subroutine ab_time_scheme_compute_coeffs
117
118end module ab_time_scheme
Interface for setting the scheme coefficients.
Adam-Bashforth scheme for time integration.
subroutine ab_time_scheme_compute_coeffs(coeffs, dt, order)
Compute the scheme coefficients.
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
Build configurations.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Base class for time integration schemes.
Utilities.
Definition utils.f90:35
Explicit Adam-Bashforth scheme for time integration.
Base abstract class for time integration schemes.