Neko  0.8.1
A portable framework for high-order spectral element flow simulations
precon_fctry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34  use precon
35  use identity
36  use device_identity
37  use jacobi
38  use sx_jacobi
39  use device_jacobi
40  use hsmg
41  use utils
42  use neko_config
43  implicit none
44 
45 contains
46 
48  subroutine precon_factory(pc, pctype)
49  class(pc_t), target, allocatable, intent(inout) :: pc
50  character(len=*), intent(in) :: pctype
51 
52  if (allocated(pc)) then
53  call precon_destroy(pc)
54  deallocate(pc)
55  end if
56 
57  if (trim(pctype) .eq. 'jacobi') then
58  if (neko_bcknd_sx .eq. 1) then
59  allocate(sx_jacobi_t::pc)
60  else if (neko_bcknd_device .eq. 1) then
61  allocate(device_jacobi_t::pc)
62  else
63  allocate(jacobi_t::pc)
64  end if
65  else if (pctype(1:4) .eq. 'hsmg') then
66  allocate(hsmg_t::pc)
67  else if(trim(pctype) .eq. 'ident') then
68  if (neko_bcknd_device .eq. 1) then
69  allocate(device_ident_t::pc)
70  else
71  allocate(ident_t::pc)
72  end if
73  else
74  call neko_error('Unknown preconditioner '//trim(pctype))
75  end if
76 
77  end subroutine precon_factory
78 
80  subroutine precon_destroy(pc)
81  class(pc_t), allocatable, intent(inout) :: pc
82 
83  if (allocated(pc)) then
84  select type(pcp => pc)
85  type is(jacobi_t)
86  call pcp%free()
87  type is(sx_jacobi_t)
88  call pcp%free()
89  type is(device_jacobi_t)
90  call pcp%free()
91  type is (hsmg_t)
92  call pcp%free()
93  end select
94  end if
95 
96  end subroutine precon_destroy
97 
98 end module precon_fctry
Identity Krylov preconditioner for accelerators.
Jacobi preconditioner accelerator backend.
Krylov preconditioner.
Definition: pc_hsmg.f90:61
Krylov preconditioner (identity)
Definition: pc_identity.f90:34
Jacobi preconditioner.
Definition: pc_jacobi.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_sx
Definition: neko_config.f90:39
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
subroutine precon_destroy(pc)
Destroy a preconditioner.
subroutine precon_factory(pc, pctype)
Create a preconditioner.
Krylov preconditioner.
Definition: precon.f90:34
Jacobi preconditioner SX-Aurora backend.
Utilities.
Definition: utils.f90:35
Defines a canonical Krylov preconditioner for accelerators.
Defines a jacobi preconditioner.
Defines a canonical Krylov preconditioner.
Definition: pc_identity.f90:42
Defines a jacobi preconditioner.
Definition: pc_jacobi.f90:45
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40
Defines a jacobi preconditioner for SX-Aurora.