Choral
cell_mod.f90
Go to the documentation of this file.
1 
80 
81 
82 module cell_mod
83 
85  use real_type
86  use basic_tools
87 
88  implicit none
89  private
90 
91  ! %----------------------------------------%
92  ! | |
93  ! | PUBLIC DATA |
94  ! | |
95  ! %----------------------------------------%
96 
97  public :: cell_init ! tested
98  !!
99  public :: cell_nbnodes, cell_nbvtx
100  public :: cell_dim , cell_name
101  public :: cell_nbfc , cell_fc_nbvtx
102  public :: cell_fc_vtx , cell_geo
103  public :: cell_nbed , cell_ed_vtx
104  public :: cell_to_gmsh, gmsh_to_cell
105  public :: cell_coord , cell_nbitf
106 
107 
110 
112  character(LEN=4), dimension(CELL_TOT_NB) :: cell_name
113 
115  integer, dimension(CELL_TOT_NB) :: cell_nbnodes=0
116 
118  integer, dimension(CELL_TOT_NB) :: cell_nbvtx =0
119 
121  integer, dimension(CELL_TOT_NB) :: cell_nbed =0
122 
124  integer, dimension(CELL_TOT_NB) :: cell_nbfc =0
125 
127  integer, dimension(CELL_TOT_NB) :: cell_nbitf =0
128 
130  integer, dimension(CELL_TOT_NB) :: cell_dim =0
131 
133  integer, dimension(CELL_TOT_NB) :: cell_geo =0
134 
137  integer, dimension(2,CELL_MAX_NBED,CELL_TOT_NB)::cell_ed_vtx=0
138 
141  integer, dimension(CELL_MAX_NBFC,CELL_TOT_NB) :: cell_fc_nbvtx =0
142 
143 
146  integer, dimension(CELL_MAX_NBVTX,CELL_MAX_NBFC,CELL_TOT_NB)::cell_fc_vtx=0
147 
149  integer, dimension(CELL_TOT_NB) :: cell_to_gmsh = -1
150 
152  integer, dimension(:), allocatable :: gmsh_to_cell
153 
155  type(r_2d), dimension(CELL_TOT_NB) :: cell_coord
156 
157 
158 contains
159 
160 
162  subroutine cell_init(b)
163  logical, intent(in) :: b
164 
165  if (b) write(*,*)"cell_mod : cell_init"
166  call def_cell()
167  call def_cell_fc()
168  call def_cell_ed_vtx()
169  call def_gmsh_to_cell()
170  call def_cell_coord()
171 
172  end subroutine cell_init
173 
177  subroutine def_cell()
179  cell_name(cell_vtx ) = 'VTX'
180  cell_name(cell_edg ) = 'EDG'
181  cell_name(cell_edg_2) = 'EDG2'
182  cell_name(cell_trg ) = 'TRG'
183  cell_name(cell_trg_2) = 'TRG2'
184  cell_name(cell_tet ) = 'TET'
185  cell_name(cell_tet_2) = 'TET2'
186 
187  cell_nbnodes(cell_vtx ) = 1
188  cell_nbnodes(cell_edg ) = 2
190  cell_nbnodes(cell_trg ) = 3
192  cell_nbnodes(cell_tet ) = 4
194 
195  if ( cell_max_nbnodes /= maxval(cell_nbnodes) ) &
196  & call quit( 'cell_mod: def_CELL: 1' )
197 
198  cell_nbvtx(cell_vtx ) = 1
199  cell_nbvtx(cell_edg ) = 2
201  cell_nbvtx(cell_trg ) = 3
203  cell_nbvtx(cell_tet ) = 4
205 
206  if (cell_max_nbvtx /= maxval(cell_nbvtx)) &
207  & call quit( 'cell_mod: def_CELL: 2' )
208 
209  cell_nbed(cell_trg ) = 3
210  cell_nbed(cell_trg_2) = 3
211  cell_nbed(cell_tet ) = 6
212  cell_nbed(cell_tet_2) = 6
213 
214  if (cell_max_nbed /= maxval(cell_nbed)) &
215  & call quit( 'cell_mod: def_CELL: 3' )
216 
217  cell_nbfc(cell_tet ) = 4
218  cell_nbfc(cell_tet_2) = 4
219  if (cell_max_nbfc /= maxval(cell_nbfc)) &
220  & call quit( 'cell_mod: def_CELL: 4' )
221 
222  cell_nbitf(cell_vtx ) = 0
223  cell_nbitf(cell_edg ) = 2
225  cell_nbitf(cell_trg ) = 3
227  cell_nbitf(cell_tet ) = 4
229  if (cell_max_nbitf /= maxval(cell_nbitf)) &
230  & call quit( 'cell_mod: def_CELL: 5' )
231 
232 
233  cell_dim(cell_vtx ) = 0
234  cell_dim(cell_edg ) = 1
235  cell_dim(cell_edg_2) = 1
236  cell_dim(cell_trg ) = 2
237  cell_dim(cell_trg_2) = 2
238  cell_dim(cell_tet ) = 3
239  cell_dim(cell_tet_2) = 3
240 
248 
249  end subroutine def_cell
250 
252  subroutine def_cell_ed_vtx()
254  integer :: c
255 
256  c = cell_trg
257  cell_ed_vtx(:,1,c) = (/1,2/)
258  cell_ed_vtx(:,2,c) = (/2,3/)
259  cell_ed_vtx(:,3,c) = (/3,1/)
260 
262 
263  c = cell_tet
264  cell_ed_vtx(:,1,c) = (/1,2/)
265  cell_ed_vtx(:,2,c) = (/2,3/)
266  cell_ed_vtx(:,3,c) = (/3,1/)
267  cell_ed_vtx(:,4,c) = (/1,4/)
268  cell_ed_vtx(:,5,c) = (/2,4/)
269  cell_ed_vtx(:,6,c) = (/3,4/)
270 
271  c = cell_tet_2
272  cell_ed_vtx(:,1,c) = (/1,2/)
273  cell_ed_vtx(:,2,c) = (/2,3/)
274  cell_ed_vtx(:,3,c) = (/3,1/)
275  cell_ed_vtx(:,4,c) = (/1,4/)
276  cell_ed_vtx(:,5,c) = (/2,4/)
277  cell_ed_vtx(:,6,c) = (/3,4/)
278 
279  end subroutine def_cell_ed_vtx
280 
281 
283  subroutine def_cell_fc()
285  integer :: c
286 
287  c = cell_tet
288  cell_fc_nbvtx(1:4, c) = 3
289 
290  c = cell_tet_2
291  cell_fc_nbvtx(1:4, c) = 3
292 
293  if (cell_max_fc_nbvtx /= maxval(cell_fc_nbvtx)) &
294  & call quit( 'cell_mod: def_CELL_FC' )
295 
296  c = cell_tet
297  cell_fc_vtx(1:3,1,c) = (/1,2,3/)
298  cell_fc_vtx(1:3,2,c) = (/1,2,4/)
299  cell_fc_vtx(1:3,3,c) = (/2,3,4/)
300  cell_fc_vtx(1:3,4,c) = (/3,1,4/)
301 
302  c = cell_tet_2
303  cell_fc_vtx(1:3,1,c) = (/1,2,3/)
304  cell_fc_vtx(1:3,2,c) = (/1,2,4/)
305  cell_fc_vtx(1:3,3,c) = (/2,3,4/)
306  cell_fc_vtx(1:3,4,c) = (/3,1,4/)
307 
308  end subroutine def_cell_fc
309 
311  subroutine def_gmsh_to_cell()
313  integer :: ii, jj
314 
315  cell_to_gmsh(cell_vtx) = 15
322 
323  call allocmem( gmsh_to_cell, maxval(cell_to_gmsh) )
324 
325  gmsh_to_cell = -1
326  do ii=1, cell_tot_nb
327 
328  jj = cell_to_gmsh(ii)
329  gmsh_to_cell(jj) = ii
330 
331  end do
332 
333  end subroutine def_gmsh_to_cell
334 
335 
337  subroutine def_cell_coord()
339  integer :: ct
340 
341  do ct=1, cell_tot_nb
342  cell_coord(ct) = r_2d( cell_dim(ct), cell_nbnodes(ct) )
343  cell_coord(ct)%y = 0.0_rp
344  end do
345 
346  ct = cell_edg
347  cell_coord(ct)%y(1, 2) = 1.0_rp
348 
349  ct = cell_trg
350  cell_coord(ct)%y(1, 2) = 1.0_rp
351  cell_coord(ct)%y(2, 3) = 1.0_rp
352 
353  ct = cell_tet
354  cell_coord(ct)%y(1, 2) = 1.0_rp
355  cell_coord(ct)%y(2, 3) = 1.0_rp
356  cell_coord(ct)%y(3, 4) = 1.0_rp
357 
358  ct = cell_tet_2
359  cell_coord(ct)%y(1, 2) = 1.0_rp
360  cell_coord(ct)%y(2, 3) = 1.0_rp
361  cell_coord(ct)%y(3, 4) = 1.0_rp
362 
363  cell_coord(ct)%y( 1, 5) = 0.5_rp
364  cell_coord(ct)%y(1:2, 6) = 0.5_rp
365  cell_coord(ct)%y( 2, 7) = 0.5_rp
366 
367  cell_coord(ct)%y( 3, 8:10) = 0.5_rp
368  cell_coord(ct)%y( 2, 9 ) = 0.5_rp
369  cell_coord(ct)%y( 1, 10 ) = 0.5_rp
370 
371  ct = cell_edg_2
372  cell_coord(ct)%y(:, 1:2) = cell_coord(cell_edg)%y
373  cell_coord(ct)%y(1, 3 ) = 0.5_rp
374 
375  ct = cell_trg_2
376  cell_coord(ct)%y(:, 1:3) = cell_coord(cell_trg)%y
377  cell_coord(ct)%y(:, 4 ) = (/ 0.5_rp, 0.0_rp /)
378  cell_coord(ct)%y(:, 5 ) = (/ 0.5_rp, 0.5_rp /)
379  cell_coord(ct)%y(:, 6 ) = (/ 0.0_rp, 0.5_rp /)
380 
381  end subroutine def_cell_coord
382 
383 
384 end module cell_mod
integer, parameter, public cell_max_nbed
Cell maximal number of edges.
integer, dimension(:), allocatable, public gmsh_to_cell
Dictionnary GMSH cell desc. –> CHORAL cell desc.
Definition: cell_mod.f90:152
integer, parameter cell_edg_2
Quadratic edge.
integer, parameter cell_tot_nb
Number of CELL types.
integer, parameter cell_tet
Tetrahedron.
BASIC TOOLS
Definition: basic_tools.F90:8
integer, dimension(cell_tot_nb), public cell_nbvtx
Number of vertexes for each cell type.
Definition: cell_mod.f90:118
integer, dimension(cell_tot_nb), public cell_nbfc
Number of faces for each cell type.
Definition: cell_mod.f90:124
subroutine, public cell_init(b)
Initialisation of all the arrays CELL_XXX.
Definition: cell_mod.f90:163
subroutine, public quit(message)
Stop program execution, display an error messahe.
integer, dimension(cell_tot_nb), public cell_geo
Associated reference cell for each cell type.
Definition: cell_mod.f90:133
integer, parameter cell_edg
Edge (line segment)
integer, dimension(cell_tot_nb), public cell_dim
Dimension for each cell type.
Definition: cell_mod.f90:130
integer, dimension(cell_max_nbvtx, cell_max_nbfc, cell_tot_nb), public cell_fc_vtx
CELL_FC_VTX(1:n, fc, cl) = vertexes for the face fc of the cell cl, with n = CELL_FC_NBVTX(fc, cl)
Definition: cell_mod.f90:146
REAL NUMBERS PRECISION IN CHORAL: selects simple/double/quad
Definition: real_type.F90:33
integer, dimension(cell_max_nbfc, cell_tot_nb), public cell_fc_nbvtx
CELL_FC_NBVTX(fc, cl) = number of vertexes for the face fc of the cell cl.
Definition: cell_mod.f90:141
integer, dimension(cell_tot_nb), public cell_nbnodes
Number of nodes for each cell type.
Definition: cell_mod.f90:115
type(r_2d), dimension(cell_tot_nb), public cell_coord
Cell node coordinates.
Definition: cell_mod.f90:155
integer, parameter cell_trg
Triangle.
integer, dimension(cell_tot_nb), public cell_nbitf
Number of interfaces for each cell type.
Definition: cell_mod.f90:127
integer, parameter cell_tet_2
Quadratic tetrahedron.
CHORAL CONSTANTS
integer, parameter, public cell_max_fc_nbvtx
Face maximal number of vertexes.
integer, parameter cell_trg_2
Quadratic triangle.
integer, dimension(cell_tot_nb), public cell_nbed
Number of edges for each cell type.
Definition: cell_mod.f90:121
integer, parameter, public cell_max_nbfc
Cell maximal number of faces.
integer, parameter, public cell_max_nbitf
Cell maximal number of interfaces.
allocate memory for real(RP) arrays
Definition: real_type.F90:158
subroutine def_cell_coord()
Array CELL_COORD initialisation.
Definition: cell_mod.f90:338
subroutine def_cell_ed_vtx()
Array CELL_ED_VTX initialisation.
Definition: cell_mod.f90:253
character(len=4), dimension(cell_tot_nb), public cell_name
CELL_ARRAY Arrays describing cells.
Definition: cell_mod.f90:112
integer, dimension(cell_tot_nb), public cell_to_gmsh
Dictionnary cell desc. –> GMSH cell desc.
Definition: cell_mod.f90:149
integer, parameter cell_vtx
Vertex.
R_2D: this type will allow to define arrays of 2D real arrays.
Definition: real_type.F90:114
subroutine def_cell()
Initialisation of : CELL_NAME, CELL_NBNODES, CELL_NBVTX, CELL_NBED, CELL_NBFC , CELL_NBITF, CELL_DIM , CELL_GEO.
Definition: cell_mod.f90:178
integer, dimension(2, cell_max_nbed, cell_tot_nb), public cell_ed_vtx
CELL_ED_VTX(1:2, ed, cl) = vertexes for the edge ed of the cell cl.
Definition: cell_mod.f90:137
subroutine def_cell_fc()
Array CELL_FC initialisation.
Definition: cell_mod.f90:284
DEFINITION OF GEOMETRICAL CELLS (for meshes)
Definition: cell_mod.f90:82
subroutine def_gmsh_to_cell()
Array CELL_TO_GMSH, GMSH_TO_CELL initialisation.
Definition: cell_mod.f90:312
integer, parameter, public cell_max_nbvtx
Cell maximal number of vertexes.
integer, parameter, public cell_max_nbnodes
Cell maximal number of nodes.