Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
navier6.F90
Go to the documentation of this file.
1 !===============================================================================
2 ! pff@cfm.brown.edu 3/19/96
3 
4 
6 ! problems with finite elements on distributed memory machines.
7 
8 ! The overall hierarchy of data structures is as follows:
9 
10 ! System - index set denoted by _glob
11 
12 ! Processor - index set denoted by _pglob
13 
14 ! .Domain - index set denoted by _dloc (1,2,3,...,n_dloc)
15 
16 ! .Sp.Elem - index set denoted by _sloc (1,2,3,...,n_sloc)
17 
18 
19 ! A critical component of the parallel DD solver is the gather-scatter
20 ! operation. As there may be more than one domain on a processor, and
21 ! communication must be minimized, it is critical that communication be
22 ! processor oriented, not domain oriented. Hence domains will access data
23 ! via the dloc_to_pglob interface, and the pglob indexed variables will
24 ! be accessed via a gather-scatter which interacts via the pglob_glob
25 ! interface. Noticed that, in a uni-processor application, the pglob_glob
26 ! map will be simply the identity.
27 
28 !===============================================================================
29 
31 subroutine set_overlap
32  use size_m, only : lx1, nid
33  use input, only : param, ifaxis, ifmgrid, ifmhd, ifsplit, ifldmhd
34  use tstep, only : ifield
35  use hsmg, only : use_spectral_coarse
37  implicit none
38 
39 !max real(DP) :: x(2*ltotd), y(2*ltotd), z(2*ltotd)
40 
41  integer :: npass, ipass
42 
43  if (lx1 == 2) param(43)=1.
44  if (lx1 == 2 .AND. nid == 0) write(6,*) 'No mgrid for lx1=2!'
45 
46  if (ifaxis) ifmgrid = .false.
47  if (param(43) /= 0) ifmgrid = .false.
48 
49  npass = 1
50  if (ifmhd) npass = 2
51  do ipass=1,npass
52  ifield = 1
53 
54  if (ifsplit .AND. ifmgrid) then
55  if (ipass > 1) ifield = ifldmhd
56 
57  call swap_lengths
58  call gen_fast_spacing()
59 
60  call hsmg_setup
61  call h1mg_setup
62 
63  elseif ( .NOT. ifsplit) then ! Pn-Pn-2
64 #if 0
65  if (ipass > 1) ifield = ifldmhd
66 
67  if (param(44) == 1) then ! Set up local overlapping solves
68  call set_fem_data_l2(nel_proc,ndom,n_o,x,y,z,tri)
69  else
70  call swap_lengths
71  endif
72 
73  e = 1
74  if (ifield > 1) e = nelv+1
75 
76  call gen_fast_spacing(x,y,z)
77  call gen_fast(df(1,e),sr(1,e),ss(1,e),st(1,e),x,y,z)
78 
79  call init_weight_op
80  if (param(43) == 0) call hsmg_setup
81 #endif
82  endif
83 
84  if (.not. use_spectral_coarse) call set_up_h1_crs
85 
86  enddo
87 
88  return
89 end subroutine set_overlap
90 !-----------------------------------------------------------------------
subroutine, public hsmg_setup()
Sets up hybrid Schwarz multi-grid preconditioner.
Definition: hsmg.F90:50
cleaned
Definition: tstep_mod.F90:2
Input parameters from preprocessors.
Definition: input_mod.F90:11
subroutine gen_fast_spacing()
Generate fast diagonalization matrices for each element.
Definition: fast3d.F90:3
Module containing data for HSMG.
Definition: hsmg_mod.F90:3
subroutine swap_lengths()
Definition: fast3d.F90:453
subroutine, public h1mg_setup()
Sets up Poisson preconditioner.
Definition: hsmg.F90:74
static double y[NR *NS *NT *N]
Definition: obbox_test.c:31
static double z[NR *NS *NT *N]
Definition: obbox_test.c:31