FE-Project
Loading...
Searching...
No Matches
scale_mesh_linedom1d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9!-------------------------------------------------------------------------------
10#include "scaleFElib.h"
12
13 !-----------------------------------------------------------------------------
14 !
15 !++ used modules
16 !
17 use scale_precision
18
19 use scale_mesh_base1d, only: &
22
26
27 !-----------------------------------------------------------------------------
28 implicit none
29 private
30
31 !-----------------------------------------------------------------------------
32 !
33 !++ Public type & procedure
34 !
35 type, extends(meshbase1d), public :: meshlinedom1d
36 contains
37 procedure :: init => meshlinedom1d_init
38 procedure :: final => meshlinedom1d_final
39 procedure :: generate => meshlinedom1d_generate
40 end type meshlinedom1d
41
44
45 !-----------------------------------------------------------------------------
46 !
47 !++ Public parameters & variables
48 !
49
50 !-----------------------------------------------------------------------------
51 !
52 !++ Private procedure
53 !
54
55 !-----------------------------------------------------------------------------
56 !
57 !++ Private parameters & variables
58 !
59
60contains
61 subroutine meshlinedom1d_init(this, & ! (inout)
62 neg, & ! (in)
63 dom_xmin, dom_xmax, & ! (in)
64 refelem, nlocalmeshperprc, & ! (in)
65 nproc, myrank, fx ) ! (in)
66
67 implicit none
68
69 class(meshlinedom1d), intent(inout) :: this
70 integer, intent(in) :: neg
71 real(rp), intent(in) :: dom_xmin
72 real(rp), intent(in) :: dom_xmax
73 type(lineelement), intent(in), target :: refelem
74 integer, intent(in) :: nlocalmeshperprc
75 integer, intent(in), optional :: nproc
76 integer, intent(in), optional :: myrank
77 real(rp), intent(in), optional :: fx(neg+1)
78 !-----------------------------------------------------------------------------
79
80 call meshbase1d_init(this, &
81 neg, &
82 dom_xmin, dom_xmax, &
83 refelem, nlocalmeshperprc, &
84 nproc, myrank, fx )
85
86 this%dom_vol = (this%xmax_gl - this%xmin_gl)
87
88 end subroutine meshlinedom1d_init
89
90 subroutine meshlinedom1d_final( this ) ! (inout)
91 implicit none
92 class(meshlinedom1d), intent(inout) :: this
93 !-----------------------------------------------------------------------------
94
95 call meshbase1d_final(this)
96
97 end subroutine meshlinedom1d_final
98
99 subroutine meshlinedom1d_generate( this ) ! (inout)
100 implicit none
101 class(meshlinedom1d), intent(inout), target :: this
102
103 integer :: n
104 integer :: p
105 type(localmesh1d), pointer :: mesh
106
107 integer :: tileid_table(this%local_mesh_num, this%prc_num)
108 integer :: panelid_table(this%local_mesh_num*this%prc_num)
109 integer :: pi_table(this%local_mesh_num*this%prc_num)
110
111 !integer :: TILE_NUM_PER_PANEL
112 integer :: tileid
113
114 !-----------------------------------------------------------------------------
115
116 ! TILE_NUM_PER_PANEL = this%LOCAL_MESH_NUM_global / 1
117
118
119 !--- Construct the connectivity of patches (only master node)
120
121 call meshbase1d_assigndomid( this, & ! (in)
122 tileid_table, panelid_table, & ! (out)
123 pi_table ) ! (out)
124
125 !--- Setup local meshes managed by my process
126
127 do n=1, this%LOCAL_MESH_NUM
128 mesh => this%lcmesh_list(n)
129 tileid = tileid_table(n, mesh%PRC_myrank+1)
130 call meshbase1d_setuplocaldom( mesh, & ! (inout)
131 tileid, panelid_table(tileid), & ! (in)
132 pi_table(tileid), this%Nprc, & ! (in)
133 this%xmin_gl, this%xmax_gl, & ! (in)
134 this%NeG / this%Nprc, this%FX(:) ) ! (in)
135
136 !---
137 ! write(*,*) "** my_rank=", mesh%PRC_myrank
138 ! write(*,*) " tileID:", mesh%tileID
139 ! write(*,*) " pnlID:", mesh%panelID, "-- i (within a panel)=", pi_table(tileID)
140 ! write(*,*) " local mesh:", n, "( total", this%LOCAL_MESH_NUM, ")"
141 ! write(*,*) " panel_connect:", this%tilePanelID_globalMap(:,mesh%tileID)
142 ! write(*,*) " tile_connect:", this%tileID_globalMap(:,mesh%tileID)
143 ! write(*,*) " face_connect:", this%tileFaceID_globalMap(:,mesh%tileID)
144 ! write(*,*) " domain size"
145 ! write(*,*) " NeX:", mesh%Ne
146 ! write(*,*) " [X]:", mesh%xmin, mesh%xmax
147 end do
148
149 this%isGenerated = .true.
150
151 return
152 end subroutine meshlinedom1d_generate
153
154 !-------------------------------------------------
155end module scale_mesh_linedom1d
module FElib / Element / Base
module FElib / Element / line
module FElib / Mesh / Local 1D
module FElib / Mesh / Base 1D
subroutine, public meshbase1d_init(this, neg, dom_xmin, dom_xmax, refelem, nlocalmeshperprc, nprocs, myrank, fx)
subroutine, public meshbase1d_final(this)
subroutine, public meshbase1d_assigndomid(this, tileid_table, panelid_table, pi_table)
subroutine, public meshbase1d_setgeometricinfo(lcmesh)
subroutine, public meshbase1d_setuplocaldom(lcmesh, tileid, panelid, i, nprc, dom_xmin, dom_xmax, ne, fx)
module FElib / Mesh / 1D domain
subroutine, public meshlinedom1d_generate(this)
subroutine, public meshlinedom1d_final(this)
subroutine, public meshlinedom1d_init(this, neg, dom_xmin, dom_xmax, refelem, nlocalmeshperprc, nproc, myrank, fx)