Choral
quadMesh_mod.f90
Go to the documentation of this file.
1 !! \li Composite integration domain:
109 
111 
112  use choral_constants
113  use choral_variables
114  use abstract_interfaces, only: r3tor
115  use basic_tools
116  use cell_mod
117  use quad_mod
118  use mesh_mod
119 
120  implicit none
121  private
122 
123 
124  ! %----------------------------------------%
125  ! | |
126  ! | PUBLIC DATA |
127  ! | |
128  ! %----------------------------------------%
129  public :: quadmesh
130 
131  public :: clear
132  public :: set
133  public :: assemble
134  public :: print
135  public :: valid
136 
137  ! %----------------------------------------%
138  ! | |
139  ! | DERIVED TYPE |
140  ! | |
141  ! %----------------------------------------%
146  type quadmesh
147 
149  type(mesh), pointer :: mesh=>null()
150 
152  integer, dimension(:), allocatable :: quadtype
153 
156  integer, dimension(0:QUAD_TOT_NB) :: quad_count = 0
157 
158  contains
159 
161  final :: quadmesh_clear
162 
163  end type quadmesh
164 
165 
166  ! %----------------------------------------%
167  ! | |
168  ! | GENERIc SUBROUTINES |
169  ! | |
170  ! %----------------------------------------%
171  interface clear
172  module procedure quadmesh_clear
173  end interface clear
174 
175  interface quadmesh
176  module procedure quadmesh_create
177  end interface quadmesh
178 
181  interface set
182  module procedure quadmesh_set
183  end interface set
184 
185  interface assemble
186  module procedure quadmesh_assemble
187  end interface assemble
188 
189  interface print
190  module procedure quadmesh_print
191  end interface print
192 
193  interface valid
194  module procedure quadmesh_valid
195  end interface valid
196 
197 contains
198 
200  subroutine quadmesh_clear(qdm)
201  type(quadmesh), intent(inout) :: qdm
202 
203  qdm%mesh => null()
204  call freemem(qdm%quadType)
205 
206  qdm%quad_count = 0
207 
208  end subroutine quadmesh_clear
209 
219  function quadmesh_create(m) result(qdm)
220  type(quadmesh) :: qdm
221  type(mesh) , target :: m
222 
223  if (choral_verb>0) write(*,*) &
224  &'quadMesh_mod : create'
225 
226  if(.NOT.valid(m)) call quit("quadMesh_mod:&
227  & quadMesh_create: mesh 'm' not valid" )
228 
229  call clear(qdm)
230  qdm%mesh => m
231  call allocmem(qdm%quadType, m%nbCl)
232  qdm%quadType = quad_none
233  if (choral_verb>1) write(*,*) &
234  &' Number of cells =', qdm%mesh%nbCl
235 
236  end function quadmesh_create
237 
238 
239 
264  subroutine quadmesh_set(qdm, quadType, f)
265  type(quadmesh), intent(inout) :: qdm
266  integer , intent(in) :: quadType
267  procedure(r3tor) , optional :: f
268 
269  logical, dimension(:), allocatable :: flag
270  integer :: cpt
271 
272  if (choral_verb>0) write(*,*) 'quadMesh_mod : quadMesh_set'
273 
274  if ( (quadtype<0) .OR. ( quadtype > quad_tot_nb) ) &
275  & call quit('quadMesh_mod: quadMesh_set:&
276  & wrong quadrature type')
277 
278  call flag_mesh_cells(flag, qdm%mesh, quad_dim(quadtype), f)
279 
280  call set_flagged_cells(qdm, cpt, quadtype, flag)
281 
282  if (choral_verb>0) write(*,*) ' ', quad_name(quadtype),&
283  & ' =',&
284  & cpt, ' cells'
285 
286  end subroutine quadmesh_set
287 
288 
289 
297  subroutine set_flagged_cells(qdm, cpt, quadType, flag)
298  type(quadmesh) , intent(inout) :: qdm
299  integer , intent(out) :: cpt
300  integer , intent(in) :: quadType
301  logical, dimension(:), intent(in) :: flag
302 
303  integer :: ii, jj, geo, cl_geo
304 
305  if(size(flag,1)/=qdm%mesh%nbCl) call quit("quadMesh_mod:&
306  & set_flagged cells: argument 'flag' incorrect" )
307 
308  geo = quad_geo(quadtype)
309  cpt = 0
310 
311  do ii=1, qdm%mesh%nbCl
312 
313  if (.NOT.flag(ii)) cycle
314 
315  jj = qdm%mesh%clType(ii)
316  cl_geo = cell_geo(jj)
317 
318  if (geo==cl_geo) then
319 
320  cpt = cpt + 1
321  qdm%quadType(ii) = quadtype
322 
323  end if
324  end do
325 
326  end subroutine set_flagged_cells
327 
328 
329 
332  subroutine quadmesh_assemble(qdm)
333  type(quadmesh), intent(inout) :: qdm
334 
335  integer :: cl, qt
336 
337  if (choral_verb>0) write(*,*)&
338  &"quadMesh_mod : assemble"
339 
340 
341  !! Set m%cell_count and m%dim
342  !!
343  qdm%quad_count = 0
344  do cl=1, qdm%mesh%nbCl
345 
346  qt = qdm%quadType(cl)
347  qdm%quad_count(qt) = qdm%quad_count(qt) + 1
348 
349  end do
350 
351  if (.NOT.valid(qdm)) call quit(&
352  & 'quadMesh_mod: mesh_ass_boundaryemble: assembling not valid')
353 
354  end subroutine quadmesh_assemble
355 
356 
357 
360  subroutine quadmesh_print(qdm)
361  type(quadmesh), intent(in) :: qdm
362 
363  integer :: ii
364 
365  write(*,*)"quadMesh_mod : print"
366  write(*,*) ' Number of cells per quad method'
367  do ii=0, quad_tot_nb
368  if ( qdm%quad_count(ii)>0) then
369 
370  write(*,*)" "&
371  &//quad_name(ii)&
372  &//" =", qdm%quad_count(ii)
373 
374  end if
375  end do
376 
377  end subroutine quadmesh_print
378 
381  function quadmesh_valid(qdm) result(b)
382  type(quadmesh), intent(in) :: qdm
383  logical :: b
384 
385  b = .false.
386  if (.NOT. ( associated(qdm%mesh) )) return
387  if (.NOT. valid(qdm%mesh)) return
388  if (.NOT. ( allocated(qdm%quadType) )) return
389 
390  b = ( size(qdm%quadType,1) == qdm%mesh%nbCl )
391  b = b .AND. ( size(qdm%quadType,1) == sum(qdm%quad_count) )
392 
393  end function quadmesh_valid
394 
395 
396 end module quadmesh_mod
Derived type mesh.
Definition: mesh_mod.F90:92
DEFINITION OF ABSTRACT_INTERFACES FOR THE LIBRARY CHORAL
BASIC TOOLS
Definition: basic_tools.F90:8
deallocate memory for real(RP) arrays
Definition: real_type.F90:163
logical function quadmesh_valid(qdm)
Check if assembling is correct.
QUADRATURE RULES ON REFERENCE CELLS
Definition: quad_mod.f90:31
subroutine, public quit(message)
Stop program execution, display an error messahe.
subroutine set_flagged_cells(qdm, cpt, quadType, flag)
Set the integration method &#39;qdm&#39; to &#39;quadType&#39; on all cells &#39;cl&#39; of compatible geometry that moreover...
integer, dimension(cell_tot_nb), public cell_geo
Associated reference cell for each cell type.
Definition: cell_mod.f90:133
subroutine, public flag_mesh_cells(flag, m, dim, f)
OUTPUT: flag(:) array of logical with size mnbCl.
Definition: mesh_mod.F90:622
integer, parameter quad_none
integer, dimension(quad_tot_nb), public quad_dim
Reference cell dimension for each quad method.
Definition: quad_mod.f90:79
DERIVED TYPE quadMesh: integration methods on meshes
integer, parameter quad_tot_nb
Total number of quadrature rules.
subroutine quadmesh_clear(qdm)
Destructor for quadMesh type.
CHORAL CONSTANTS
type(quadmesh) function quadmesh_create(m)
Constructor for the type quadMesh
DERIVED TYPE mesh
Definition: mesh_mod.F90:57
allocate memory for real(RP) arrays
Definition: real_type.F90:158
integer choral_verb
Verbosity level.
Sets an integration method quadMesh.
DEFINITION OF GLOBAL VARIABLES FOR THE LIBRARY CHORAL
subroutine quadmesh_assemble(qdm)
Sets mdim, qdmquad_count.
The type quadMesh defines integration methods on meshes.
character(len=13), dimension(0:quad_tot_nb), public quad_name
QUAD_ARRAY.
Definition: quad_mod.f90:70
subroutine quadmesh_print(qdm)
number of cells per fe type
DEFINITION OF GEOMETRICAL CELLS (for meshes)
Definition: cell_mod.f90:82
integer, dimension(quad_tot_nb), public quad_geo
Reference cell geometry for each quad method.
Definition: quad_mod.f90:76
subroutine quadmesh_set(qdm, quadType, f)
Sets an integration method quadMesh