FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
scale_mesh_linedom1d Module Reference

module FElib / Mesh / 1D domain More...

Data Types

type  meshlinedom1d
 

Functions/Subroutines

subroutine, public meshlinedom1d_init (this, neg, dom_xmin, dom_xmax, refelem, nlocalmeshperprc, nproc, myrank, fx)
 
subroutine, public meshlinedom1d_final (this)
 
subroutine, public meshlinedom1d_generate (this)
 

Detailed Description

module FElib / Mesh / 1D domain

Description
Manage mesh data of 1D domain for element-based methods
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshlinedom1d_init()

subroutine, public scale_mesh_linedom1d::meshlinedom1d_init ( class(meshlinedom1d), intent(inout) this,
integer, intent(in) neg,
real(rp), intent(in) dom_xmin,
real(rp), intent(in) dom_xmax,
type(lineelement), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in), optional nproc,
integer, intent(in), optional myrank,
real(rp), dimension(neg+1), intent(in), optional fx )

Definition at line 61 of file scale_mesh_linedom1d.F90.

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

References scale_mesh_base1d::meshbase1d_init().

Referenced by scale_mesh_linedom1d::meshlinedom1d::generate().

◆ meshlinedom1d_final()

subroutine, public scale_mesh_linedom1d::meshlinedom1d_final ( class(meshlinedom1d), intent(inout) this)

Definition at line 90 of file scale_mesh_linedom1d.F90.

91 implicit none
92 class(MeshLineDom1D), intent(inout) :: this
93 !-----------------------------------------------------------------------------
94
95 call meshbase1d_final(this)
96

References scale_mesh_base1d::meshbase1d_final().

Referenced by scale_mesh_linedom1d::meshlinedom1d::generate().

◆ meshlinedom1d_generate()

subroutine, public scale_mesh_linedom1d::meshlinedom1d_generate ( class(meshlinedom1d), intent(inout), target this)

Definition at line 99 of file scale_mesh_linedom1d.F90.

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

References scale_mesh_base1d::meshbase1d_assigndomid(), and scale_mesh_base1d::meshbase1d_setuplocaldom().

Referenced by scale_mesh_linedom1d::meshlinedom1d::generate().