Choral
All Classes Namespaces Files Functions Variables Modules
fe_mod.F90
Go to the documentation of this file.
1 
18 
19 module fe_mod
20 
22  use real_type
23  use basic_tools
24  use fe_func_mod
25  use cell_mod
26 
27  implicit none
28  private
29 
30 
31  ! %----------------------------------------%
32  ! | |
33  ! | PUBLIC DATA |
34  ! | |
35  ! %----------------------------------------%
36 
37  public :: fe_init ! TESTED
38  public :: scal_fe ! TESTED
39  public :: scalgrad_fe ! TESTED
40  public :: vect_fe ! TESTED
41  public :: vectdiv_fe ! TESTED
42 
43  public :: fe_geo, fe_nbdof
44  public :: fe_dim, fe_name
45  public :: fe_dof_geo , fe_dof_geoelmt
46  public :: fe_dof_nbvtx, fe_dof_vtx
47  public :: fe_dof_bary , fe_dof_type
48  public :: fe_func
49  public :: fe_dof_coord
50 
51 
52  ! %----------------------------------------%
53  ! | |
54  ! | FE_XXS ARRAY DEFINITION |
55  ! | |
56  ! %----------------------------------------%
57 
60 
62  type(fe_functions), dimension(FE_TOT_NB) :: fe_func
63 
65  character(LEN=16), dimension(0:FE_TOT_NB) :: fe_name
66 
68  integer , dimension(FE_TOT_NB):: fe_nbdof
69 
71  integer , dimension(FE_TOT_NB):: fe_geo
72 
74  integer , dimension(FE_TOT_NB):: fe_dim
75 
78  integer , dimension(FE_MAX_NBDOF,FE_TOT_NB):: fe_dof_geo = -1
79 
82  integer , dimension(FE_MAX_NBDOF,FE_TOT_NB):: fe_dof_type = -1
83 
99  integer , dimension(FE_MAX_NBDOF,FE_TOT_NB):: fe_dof_geoelmt= -1
100 
104  integer , dimension(FE_MAX_NBDOF,FE_TOT_NB):: fe_dof_nbvtx = -1
105 
112  integer ,dimension(CELL_MAX_NBVTX,FE_MAX_NBDOF,FE_TOT_NB)::fe_dof_vtx =-1
113 
121  real(RP),dimension(CELL_MAX_NBVTX,FE_MAX_NBDOF,FE_TOT_NB)::fe_dof_bary=0._rp
122 
125  type(r_2d), dimension(FE_TOT_NB) :: fe_dof_coord
126 
127 
128 
129 contains
130 
135  subroutine scal_fe(val, nbDof, x, dim, ft)
136  integer , intent(in) :: dim, nbDof, ft
137  real(RP), dimension(nbDof), intent(out) :: val
138  real(RP), dimension(dim) , intent(in) :: x
139 
140 #if DBG>0
141  if ( (ft<=0) .OR. (ft>fe_tot_nb) ) &
142  & call quit( 'fe_mod: scal_fe:&
143  & wrong fe type' )
144  if ( .NOT. associated( fe_func(ft)%u) ) &
145  & call quit( "fe_mod: scal_fe:&
146  & scalar fe function 'FE_FUNC(ft)%u' not defined " )
147  if ( fe_nbdof(ft) /= nbdof ) &
148  & call quit( "fe_mod: scal_fe:&
149  & argument 'nbDof' incorrect" )
150  if ( fe_dim(ft) /= dim ) &
151  & call quit( "fe_mod: scal_fe:&
152  & argument 'dim' incorrect")
153 #endif
154 
155  call fe_func(ft)%u( val, nbdof, x, dim)
156 
157  end subroutine scal_fe
158 
159 
165  subroutine scalgrad_fe(val, nbDof, x, dim, ft)
166  integer , intent(in) :: dim, nbDof, ft
167  real(RP), dimension(dim, nbDof), intent(out) :: val
168  real(RP), dimension(dim) , intent(in) :: x
169 
170 #if DBG>0
171  if ( (ft<=0) .OR. (ft>fe_tot_nb) ) &
172  & call quit( "fe_mod: scalGrad_fe:&
173  & wrong fe type" )
174  if ( .NOT. associated( fe_func(ft)%grad_u) ) &
175  & call quit( "fe_mod: scalGrad_fe:&
176  & gradient fe function 'FE_FUNC(ft)%grad_u' not defined" )
177  if ( fe_nbdof(ft) /= nbdof ) &
178  & call quit( "fe_mod: scalGrad_fe:&
179  & argument 'nbDof' incorrect" )
180  if ( fe_dim(ft) /= dim ) &
181  & call quit( "fe_mod: scalGrad_fe:&
182  & argument 'dim' incorrect" )
183 #endif
184 
185  call fe_func(ft)%grad_u( val, nbdof, x, dim)
186 
187  end subroutine scalgrad_fe
188 
189 
194  subroutine vect_fe(val, nbDof, x, dim, ft)
195  integer , intent(in) :: dim, nbDof, ft
196  real(RP), dimension(dim,nbDof), intent(out) :: val
197  real(RP), dimension(dim) , intent(in) :: x
198 
199 #if DBG>0
200  if ( (ft<=0) .OR. (ft>fe_tot_nb) ) &
201  & call quit( "fe_mod: vect_fe: &
202  & wrong fe type" )
203  if ( .NOT. associated( fe_func(ft)%phi) )&
204  & call quit( "fe_mod: vect_fe: &
205  & vector fe function ' FE_FUNC(ft)%phi' not defined" )
206  if ( fe_nbdof(ft) /= nbdof ) &
207  & call quit( "fe_mod: vect_fe: &
208  & argument 'nbDof' incorrect" )
209  if ( fe_dim(ft) /= dim ) &
210  & call quit( "fe_mod: vect_fe: &
211  & argument 'dim' incorrect" )
212 #endif
213 
214  call fe_func(ft)%phi( val, nbdof, x, dim)
215 
216  end subroutine vect_fe
217 
218 
224  subroutine vectdiv_fe(val, nbDof, x, dim, ft)
225  integer , intent(in) :: dim, nbDof, ft
226  real(RP), dimension(nbDof), intent(out) :: val
227  real(RP), dimension(dim) , intent(in) :: x
228 
229 #if DBG>0
230  if ( (ft<=0) .OR. (ft>fe_tot_nb) ) &
231  & call quit( "fe_mod: vectDiv_fe:&
232  & wrong fe type" )
233  if ( .NOT. associated( fe_func(ft)%div_phi) ) &
234  & call quit( "fe_mod: vectDiv_fe:&
235  & divergence fe function 'FE_FUNC(ft)%div_phi' not defined" )
236  if ( fe_nbdof(ft) /= nbdof ) &
237  & call quit( "fe_mod: vectDiv_fe:&
238  & argument 'nbDof' incorrect" )
239  if ( fe_dim(ft) /= dim ) &
240  & call quit( "fe_mod: vectDiv_fe:&
241  & argument 'dim' incorrect" )
242 #endif
243 
244  call fe_func(ft)%div_phi( val, nbdof, x, dim)
245 
246  end subroutine vectdiv_fe
247 
248 
251  subroutine fe_init(b)
252  logical, intent(in) :: b
253 
254  if (b) write(*,*) "fe_mod : fe_init"
255 
256  !! FE CONSTANTS SET BY HAND
257  call def_fe()
258  call def_fe_dof_type()
259  call def_fe_dof_geo()
260  call def_fe_dof_geoelmt()
261  call def_fe_dof_vtx()
262  call def_fe_dof_bary()
263  call def_fe_func()
264 
265  !! FE CONSTANTS AUTOMATICALLY SET
266  call def_fe_dim()
267  call def_fe_dof_nbvtx()
268  call def_fe_dof_coord()
269 
270  end subroutine fe_init
271 
272 
275  subroutine def_fe()
277  fe_name(fe_none) = 'VOID'
278  fe_name(fe_p1_1d) = 'P1_1D'
279  fe_name(fe_p2_1d) = 'P2_1D'
280  fe_name(fe_p3_1d) = 'P3_1D'
281  fe_name(fe_p0_2d) = 'P0_2D'
282  fe_name(fe_p1_2d) = 'P1_2D'
283  fe_name(fe_p2_2d) = 'P2_2D'
284  fe_name(fe_p3_2d) = 'P3_2D'
285  fe_name(fe_p1_3d) = 'P1_3D'
286  fe_name(fe_p1_2d_disc_ortho) = 'P1_2D_DISC_ORTHO'
287  fe_name(fe_rt0_2d) = 'RT0_2D'
288  fe_name(fe_rt1_2d_2) = 'RT1_2D_2'
289  fe_name(fe_p2_3d) = 'P2_3D'
290  fe_name(fe_p0_1d) = 'P0_1D'
291  fe_name(fe_rt0_1d) = 'RT0_1D'
292  fe_name(fe_p1_1d_disc_ortho) = 'P1_1D_DISC_ORTHO'
293  fe_name(fe_rt1_1d) = 'RT1_1D'
294  fe_name(fe_rt1_1d_dual) = 'RT1_1D_DUAL'
295  fe_name(fe_rt1_1d_2) = 'RT1_1D_2'
296  fe_name(fe_rt1_1d_2_dual)= 'RT1_1D_2_DUAL'
297  fe_name(fe_rt0_1d_dual) = 'RT0_1D_DUAL'
298 
299  fe_nbdof(fe_p1_1d) = 2
300  fe_nbdof(fe_p2_1d) = 3
301  fe_nbdof(fe_p3_1d) = 4
302  fe_nbdof(fe_p0_2d) = 1
303  fe_nbdof(fe_p1_2d) = 3
304  fe_nbdof(fe_p2_2d) = 6
305  fe_nbdof(fe_p3_2d) = 10
306  fe_nbdof(fe_p1_3d) = 4
308  fe_nbdof(fe_rt0_2d) = 3
309  fe_nbdof(fe_rt1_2d_2) = 8
310  fe_nbdof(fe_p2_3d) = 10
311  fe_nbdof(fe_p0_1d) = 1
312  fe_nbdof(fe_rt0_1d)= 2
314  fe_nbdof(fe_rt1_1d)= 3
319 
320  if ( fe_max_nbdof /= maxval(fe_nbdof) ) &
321  & call quit( 'fe_mod: def_FE_NBDOF' )
322 
343 
344  end subroutine def_fe
345 
346 
349  subroutine def_fe_dof_type()
351  integer :: f, n
352 
353  f = fe_p1_1d
354  n = fe_nbdof(f)
355  fe_dof_type(1:n, f ) = fe_dof_lag
356 
357  f = fe_p2_1d
358  n = fe_nbdof(f)
359  fe_dof_type(1:n, f ) = fe_dof_lag
360 
361  f = fe_p3_1d
362  n = fe_nbdof(f)
363  fe_dof_type(1:n, f ) = fe_dof_lag
364 
365  f = fe_p0_2d
366  n = fe_nbdof(f)
367  fe_dof_type(1:n, f ) = fe_dof_lag
368 
369  f = fe_p1_2d
370  n = fe_nbdof(f)
371  fe_dof_type(1:n, f ) = fe_dof_lag
372 
374  n = fe_nbdof(f)
375  fe_dof_type(1:n, f ) = fe_dof_lag
376 
377  f = fe_p2_2d
378  n = fe_nbdof(f)
379  fe_dof_type(1:n, f ) = fe_dof_lag
380 
381  f = fe_p3_2d
382  n = fe_nbdof(f)
383  fe_dof_type(1:n, f ) = fe_dof_lag
384 
385  f = fe_p1_3d
386  n = fe_nbdof(f)
387  fe_dof_type(1:n, f ) = fe_dof_lag
388 
389  f = fe_rt0_2d
390  n = fe_nbdof(f)
391  fe_dof_type(1:n, f ) = fe_dof_flx
392 
393  f = fe_rt1_2d_2
394  n = fe_nbdof(f)
395  fe_dof_type(1:6, f ) = fe_dof_nrm_trace
396  fe_dof_type( 7, f ) = fe_dof_comp1
397  fe_dof_type( n, f ) = fe_dof_comp2
398 
399  f = fe_p2_3d
400  n = fe_nbdof(f)
401  fe_dof_type(1:n, f ) = fe_dof_lag
402 
403  f = fe_p0_1d
404  n = fe_nbdof(f)
405  fe_dof_type(1:n, f ) = fe_dof_lag
406 
407  f = fe_rt0_1d
408  n = fe_nbdof(f)
409  fe_dof_type(1:n, f ) = fe_dof_flx
410 
412  n = fe_nbdof(f)
413  fe_dof_type(1:n, f ) = fe_dof_lag
414 
415  f = fe_rt1_1d
416  n = fe_nbdof(f)
417  fe_dof_type(1:2, f ) = fe_dof_flx
418  fe_dof_type( 3, f ) = fe_dof_comp1
419 
420  f = fe_rt1_1d_dual
421  n = fe_nbdof(f)
422  fe_dof_type(1:2, f ) = fe_dof_flx
423  fe_dof_type( 3, f ) = fe_dof_comp1
424 
425  f = fe_rt1_1d_2
426  n = fe_nbdof(f)
427  fe_dof_type(1:2, f ) = fe_dof_flx
428  fe_dof_type( 3, f ) = fe_dof_comp1
429 
430  f = fe_rt1_1d_2_dual
431  n = fe_nbdof(f)
432  fe_dof_type(1:2, f ) = fe_dof_flx
433  fe_dof_type( 3, f ) = fe_dof_comp1
434 
435  f = fe_rt0_1d_dual
436  n = fe_nbdof(f)
437  fe_dof_type(1:n, f ) = fe_dof_flx
438 
439  end subroutine def_fe_dof_type
440 
441 
444  subroutine def_fe_dof_geo()
446  integer :: f, n
447 
448  f = fe_p1_1d
449  n = fe_nbdof(f)
450  fe_dof_geo(1:n, f ) = cell_vtx
451 
452  f = fe_p2_1d
453  n = fe_nbdof(f)
454  fe_dof_geo(1:n, f ) = (/cell_vtx, cell_vtx, cell_edg/)
455 
456  f = fe_p3_1d
457  n = fe_nbdof(f)
458  fe_dof_geo(1:n, f ) = (/cell_vtx, cell_vtx, cell_edg, cell_edg/)
459 
460  f = fe_p0_2d
461  n = fe_nbdof(f)
462  fe_dof_geo(1:n, f ) = cell_trg
463 
464  f = fe_p1_2d
465  n = fe_nbdof(f)
466  fe_dof_geo(1:n, f ) = cell_vtx
467 
469  n = fe_nbdof(f)
470  fe_dof_geo(1:n, f ) = cell_trg
471 
472  f = fe_p2_2d
473  n = fe_nbdof(f)
474  fe_dof_geo(1:3, f ) = cell_vtx
475  fe_dof_geo(4:n, f ) = cell_edg
476 
477  f = fe_p3_2d
478  n = fe_nbdof(f)
479  fe_dof_geo(1:3, f ) = cell_vtx
480  fe_dof_geo(4:9, f ) = cell_edg
481  fe_dof_geo(n , f ) = cell_trg
482 
483  f = fe_p1_3d
484  n = fe_nbdof(f)
485  fe_dof_geo(1:n, f ) = cell_vtx
486 
487  f = fe_rt0_2d
488  n = fe_nbdof(f)
489  fe_dof_geo(1:n, f ) = cell_edg
490 
491  f = fe_rt1_2d_2
492  n = fe_nbdof(f)
493  fe_dof_geo(1:6, f ) = cell_edg
494  fe_dof_geo(7:n, f ) = cell_trg
495 
496  f = fe_p2_3d
497  n = fe_nbdof(f)
498  fe_dof_geo(1:4, f ) = cell_vtx
499  fe_dof_geo(5:n, f ) = cell_edg
500 
501  f = fe_p0_1d
502  n = fe_nbdof(f)
503  fe_dof_geo(1:n, f ) = cell_edg
504 
505  f = fe_rt0_1d
506  n = fe_nbdof(f)
507  fe_dof_geo(1:n, f ) = cell_vtx
508 
510  n = fe_nbdof(f)
511  fe_dof_geo(1:n, f ) = cell_edg
512 
513  f = fe_rt1_1d
514  n = fe_nbdof(f)
515  fe_dof_geo(1:2, f ) = cell_vtx
516  fe_dof_geo( 3, f ) = cell_edg
517 
518  f = fe_rt1_1d_dual
519  n = fe_nbdof(f)
520  fe_dof_geo(1:2, f ) = cell_vtx
521  fe_dof_geo( 3, f ) = cell_edg
522 
523  f = fe_rt1_1d_2
524  n = fe_nbdof(f)
525  fe_dof_geo(1:2, f ) = cell_vtx
526  fe_dof_geo( 3, f ) = cell_edg
527 
528  f = fe_rt1_1d_2_dual
529  n = fe_nbdof(f)
530  fe_dof_geo(1:2, f ) = cell_vtx
531  fe_dof_geo( 3, f ) = cell_edg
532 
533  f = fe_rt0_1d_dual
534  n = fe_nbdof(f)
535  fe_dof_geo(1:n, f ) = cell_vtx
536 
537  end subroutine def_fe_dof_geo
538 
539 
548  subroutine def_fe_dof_geoelmt()
550  integer :: f, n
551 
552  f = fe_p1_1d
553  n = fe_nbdof(f)
554  fe_dof_geoelmt(1:n, f ) = (/1, 2/)
555 
556  f = fe_p2_1d
557  n = fe_nbdof(f)
558  fe_dof_geoelmt(1:n, f ) = (/1, 2, 1/)
559 
560  f = fe_p3_1d
561  n = fe_nbdof(f)
562  fe_dof_geoelmt(1:n, f ) = (/1, 2, 1, 1/)
563 
564  f = fe_p0_2d
565  n = fe_nbdof(f)
566  fe_dof_geoelmt(1:n, f ) = (/1/)
567 
568  f = fe_p1_2d
569  n = fe_nbdof(f)
570  fe_dof_geoelmt(1:n, f ) = (/1,2,3/)
571 
573  n = fe_nbdof(f)
574  fe_dof_geoelmt(1:n, f ) = (/1,1,1/)
575 
576  f = fe_p2_2d
577  n = fe_nbdof(f)
578  fe_dof_geoelmt(1:n, f ) = (/1,2,3,1,2,3/)
579 
580  f = fe_p3_2d
581  n = fe_nbdof(f)
582  fe_dof_geoelmt(1:n, f ) = (/1,2,3,1,1,2,2,3,3,1/)
583 
584  f = fe_p1_3d
585  n = fe_nbdof(f)
586  fe_dof_geoelmt(1:n, f ) = (/1,2,3,4/)
587 
588  f = fe_rt0_2d
589  n = fe_nbdof(f)
590  fe_dof_geoelmt(1:n, f ) = (/1,2,3/)
591 
592  f = fe_rt1_2d_2
593  n = fe_nbdof(f)
594  fe_dof_geoelmt(1:n, f ) = (/1,1,2,2,3,3,1,1/)
595 
596  f = fe_p2_3d
597  n = fe_nbdof(f)
598  fe_dof_geoelmt(1:n, f ) = (/1,2,3,4,1,2,3,4,5,6/)
599 
600  f = fe_p0_1d
601  n = fe_nbdof(f)
602  fe_dof_geoelmt(1:n, f ) = (/1/)
603 
604  f = fe_rt0_1d
605  n = fe_nbdof(f)
606  fe_dof_geoelmt(1:n, f ) = (/1,2/)
607 
609  n = fe_nbdof(f)
610  fe_dof_geoelmt(1:n, f ) = (/1,1/)
611 
612  f = fe_rt1_1d
613  n = fe_nbdof(f)
614  fe_dof_geoelmt(1:3, f ) = (/1,2,1/)
615 
616  f = fe_rt1_1d_dual
617  n = fe_nbdof(f)
618  fe_dof_geoelmt(1:3, f ) = (/1,2,1/)
619 
620  f = fe_rt1_1d_2
621  n = fe_nbdof(f)
622  fe_dof_geoelmt(1:3, f ) = (/1,2,1/)
623 
624  f = fe_rt1_1d_2_dual
625  n = fe_nbdof(f)
626  fe_dof_geoelmt(1:3, f ) = (/1,2,1/)
627 
628  f = fe_rt0_1d_dual
629  n = fe_nbdof(f)
630  fe_dof_geoelmt(1:n, f ) = (/1,2/)
631 
632  end subroutine def_fe_dof_geoelmt
633 
634 
641  subroutine def_fe_dof_vtx()
643  integer :: f
644 
645  f = fe_p1_1d
646  fe_dof_vtx( 1, 1, f ) = 1
647  fe_dof_vtx( 1, 2, f ) = 2
648 
649  f = fe_p2_1d
650  fe_dof_vtx( 1, 1, f ) = 1
651  fe_dof_vtx( 1, 2, f ) = 2
652  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
653 
654  f = fe_p3_1d
655  fe_dof_vtx( 1, 1, f ) = 1
656  fe_dof_vtx( 1, 2, f ) = 2
657  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
658  fe_dof_vtx( 1: 2, 4, f ) = (/ 1, 2 /)
659 
660  f = fe_p0_2d
661  fe_dof_vtx( 1: 3, 1, f ) = (/ 1, 2, 3 /)
662 
663  f = fe_p1_2d
664  fe_dof_vtx( 1, 1, f ) = 1
665  fe_dof_vtx( 1, 2, f ) = 2
666  fe_dof_vtx( 1, 3, f ) = 3
667 
669  fe_dof_vtx( 1: 3, 1, f ) = (/ 1, 2, 3 /)
670  fe_dof_vtx( 1: 3, 2, f ) = (/ 1, 2, 3 /)
671  fe_dof_vtx( 1: 3, 3, f ) = (/ 1, 2, 3 /)
672 
673  f = fe_p2_2d
674  fe_dof_vtx( 1, 1, f ) = 1
675  fe_dof_vtx( 1, 2, f ) = 2
676  fe_dof_vtx( 1, 3, f ) = 3
677  fe_dof_vtx( 1: 2, 4, f ) = (/ 1, 2 /)
678  fe_dof_vtx( 1: 2, 5, f ) = (/ 2, 3 /)
679  fe_dof_vtx( 1: 2, 6, f ) = (/ 3, 1 /)
680 
681  f = fe_p3_2d
682  fe_dof_vtx( 1, 1 , f ) = 1
683  fe_dof_vtx( 1, 2 , f ) = 2
684  fe_dof_vtx( 1, 3 , f ) = 3
685  fe_dof_vtx( 1: 2, 4 , f ) = (/ 1, 2 /)
686  fe_dof_vtx( 1: 2, 5 , f ) = (/ 1, 2 /)
687  fe_dof_vtx( 1: 2, 6 , f ) = (/ 2, 3 /)
688  fe_dof_vtx( 1: 2, 7 , f ) = (/ 2, 3 /)
689  fe_dof_vtx( 1: 2, 8 , f ) = (/ 3, 1 /)
690  fe_dof_vtx( 1: 2, 9 , f ) = (/ 3, 1 /)
691  fe_dof_vtx( 1: 3, 10, f ) = (/ 1, 2, 3 /)
692 
693  f = fe_p1_3d
694  fe_dof_vtx( 1, 1, f ) = 1
695  fe_dof_vtx( 1, 2, f ) = 2
696  fe_dof_vtx( 1, 3, f ) = 3
697  fe_dof_vtx( 1, 4, f ) = 4
698 
699  f = fe_rt0_2d
700  fe_dof_vtx( 1: 2, 1, f ) = (/ 1, 2 /)
701  fe_dof_vtx( 1: 2, 2, f ) = (/ 2, 3 /)
702  fe_dof_vtx( 1: 2, 3, f ) = (/ 3, 1 /)
703 
704  f = fe_rt1_2d_2
705  fe_dof_vtx( 1: 2, 1, f ) = (/ 1, 2 /)
706  fe_dof_vtx( 1: 2, 2, f ) = (/ 1, 2 /)
707  fe_dof_vtx( 1: 2, 3, f ) = (/ 2, 3 /)
708  fe_dof_vtx( 1: 2, 4, f ) = (/ 2, 3 /)
709  fe_dof_vtx( 1: 2, 5, f ) = (/ 3, 1 /)
710  fe_dof_vtx( 1: 2, 6, f ) = (/ 3, 1 /)
711  fe_dof_vtx( 1: 3, 7, f ) = (/ 1, 2, 3 /)
712  fe_dof_vtx( 1: 3, 8, f ) = (/ 1, 2, 3 /)
713 
714  f = fe_p2_3d
715  fe_dof_vtx( 1 , 1 , f ) = 1
716  fe_dof_vtx( 1 , 2 , f ) = 2
717  fe_dof_vtx( 1 , 3 , f ) = 3
718  fe_dof_vtx( 1 , 4 , f ) = 4
719  fe_dof_vtx( 1: 2, 5 , f ) = (/ 1, 2 /)
720  fe_dof_vtx( 1: 2, 6 , f ) = (/ 2, 3 /)
721  fe_dof_vtx( 1: 2, 7 , f ) = (/ 3, 1 /)
722  fe_dof_vtx( 1: 2, 8 , f ) = (/ 1, 4 /)
723  fe_dof_vtx( 1: 2, 9 , f ) = (/ 2, 4 /)
724  fe_dof_vtx( 1: 2, 10, f ) = (/ 3, 4 /)
725 
726  f = fe_p0_1d
727  fe_dof_vtx( 1: 2, 1, f ) = (/ 1, 2 /)
728 
729  f = fe_rt0_1d
730  fe_dof_vtx( 1 , 1 , f ) = 1
731  fe_dof_vtx( 1 , 2 , f ) = 2
732 
734  fe_dof_vtx( 1: 2, 1, f ) = (/ 1, 2 /)
735  fe_dof_vtx( 1: 2, 2, f ) = (/ 1, 2 /)
736 
737  f = fe_rt1_1d
738  fe_dof_vtx( 1 , 1, f ) = 1
739  fe_dof_vtx( 1 , 2, f ) = 2
740  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
741 
742  f = fe_rt1_1d_dual
743  fe_dof_vtx( 1 , 1, f ) = 1
744  fe_dof_vtx( 1 , 2, f ) = 2
745  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
746 
747  f = fe_rt1_1d_2
748  fe_dof_vtx( 1 , 1, f ) = 1
749  fe_dof_vtx( 1 , 2, f ) = 2
750  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
751 
752  f = fe_rt1_1d_2_dual
753  fe_dof_vtx( 1 , 1, f ) = 1
754  fe_dof_vtx( 1 , 2, f ) = 2
755  fe_dof_vtx( 1: 2, 3, f ) = (/ 1, 2 /)
756 
757  f = fe_rt0_1d_dual
758  fe_dof_vtx( 1 , 1 , f ) = 1
759  fe_dof_vtx( 1 , 2 , f ) = 2
760 
761  end subroutine def_fe_dof_vtx
762 
763 
764 
775  subroutine def_fe_dof_bary()
777  integer :: f
778  real(RP) :: x1, x2
779 
780  f = fe_p1_1d
781  fe_dof_bary( 1, 1, f ) = 1.0_rp
782  fe_dof_bary( 1, 2, f ) = 1.0_rp
783 
784  f = fe_p2_1d
785  fe_dof_bary( 1, 1, f ) = 1.0_rp
786  fe_dof_bary( 1, 2, f ) = 1.0_rp
787  fe_dof_bary( 1: 2, 3, f ) = 0.5_rp
788 
789  f = fe_p3_1d
790  fe_dof_bary( 1, 1, f ) = 1.0_rp
791  fe_dof_bary( 1, 2, f ) = 1.0_rp
792  fe_dof_bary( 1: 2, 3, f ) = (/ re(2,3), re(1,3) /)
793  fe_dof_bary( 1: 2, 4, f ) = (/ re(1,3), re(2,3) /)
794 
795  f = fe_p0_2d
796  fe_dof_bary( 1: 3, 1, f ) = re(1,3)
797 
798  f = fe_p1_2d
799  fe_dof_bary( 1, 1, f ) = 1.0_rp
800  fe_dof_bary( 1, 2, f ) = 1.0_rp
801  fe_dof_bary( 1, 3, f ) = 1.0_rp
802 
804  fe_dof_bary( 1: 3, 1, f ) = (/ re(2,3), re(1,6), re(1,6) /)
805  fe_dof_bary( 1: 3, 2, f ) = (/ re(1,6), re(2,3), re(1,6) /)
806  fe_dof_bary( 1: 3, 3, f ) = (/ re(1,6), re(1,6), re(2,3) /)
807 
808  f = fe_p2_2d
809  fe_dof_bary( 1, 1, f ) = 1.0_rp
810  fe_dof_bary( 1, 2, f ) = 1.0_rp
811  fe_dof_bary( 1, 3, f ) = 1.0_rp
812  fe_dof_bary( 1: 2, 4, f ) = 0.5_rp
813  fe_dof_bary( 1: 2, 5, f ) = 0.5_rp
814  fe_dof_bary( 1: 2, 6, f ) = 0.5_rp
815 
816  f = fe_p3_2d
817  fe_dof_bary( 1, 1 , f ) = 1.0_rp
818  fe_dof_bary( 1, 2 , f ) = 1.0_rp
819  fe_dof_bary( 1, 3 , f ) = 1.0_rp
820  fe_dof_bary( 1: 2, 4 , f ) = (/ re(2,3), re(1,3) /)
821  fe_dof_bary( 1: 2, 5 , f ) = (/ re(1,3), re(2,3) /)
822  fe_dof_bary( 1: 2, 6 , f ) = (/ re(2,3), re(1,3) /)
823  fe_dof_bary( 1: 2, 7 , f ) = (/ re(1,3), re(2,3) /)
824  fe_dof_bary( 1: 2, 8 , f ) = (/ re(2,3), re(1,3) /)
825  fe_dof_bary( 1: 2, 9 , f ) = (/ re(1,3), re(2,3) /)
826  fe_dof_bary( 1: 3, 10, f ) = re(1,3)
827 
828  f = fe_p1_3d
829  fe_dof_bary( 1, 1, f ) = 1.0_rp
830  fe_dof_bary( 1, 2, f ) = 1.0_rp
831  fe_dof_bary( 1, 3, f ) = 1.0_rp
832  fe_dof_bary( 1, 4, f ) = 1.0_rp
833 
834  f = fe_rt0_2d
835  fe_dof_bary( 1: 2, 1, f ) = 0.5_rp
836  fe_dof_bary( 1: 2, 2, f ) = 0.5_rp
837  fe_dof_bary( 1: 2, 3, f ) = 0.5_rp
838 
839  f = fe_rt1_2d_2
840  x1 = 0.788675134594813_rp
841  x2 = 0.21132486540518702_rp
842  fe_dof_bary( 1: 2, 1, f ) = (/x1, x2/)
843  fe_dof_bary( 1: 2, 2, f ) = (/x2, x1/)
844  fe_dof_bary( 1: 2, 3, f ) = (/x1, x2/)
845  fe_dof_bary( 1: 2, 4, f ) = (/x2, x1/)
846  fe_dof_bary( 1: 2, 5, f ) = (/x1, x2/)
847  fe_dof_bary( 1: 2, 6, f ) = (/x2, x1/)
848  fe_dof_bary( 1: 3, 7, f ) = re(1,3)
849  fe_dof_bary( 1: 3, 8, f ) = re(1,3)
850 
851  f = fe_p2_3d
852  fe_dof_bary( 1 , 1 , f ) = 1.0_rp
853  fe_dof_bary( 1 , 2 , f ) = 1.0_rp
854  fe_dof_bary( 1 , 3 , f ) = 1.0_rp
855  fe_dof_bary( 1 , 4 , f ) = 1.0_rp
856  fe_dof_bary( 1: 2, 5 , f ) = 0.5_rp
857  fe_dof_bary( 1: 2, 6 , f ) = 0.5_rp
858  fe_dof_bary( 1: 2, 7 , f ) = 0.5_rp
859  fe_dof_bary( 1: 2, 8 , f ) = 0.5_rp
860  fe_dof_bary( 1: 2, 9 , f ) = 0.5_rp
861  fe_dof_bary( 1: 2, 10, f ) = 0.5_rp
862 
863  f = fe_p0_1d
864  fe_dof_bary( 1: 3, 1, f ) = 0.5_rp
865 
866  f = fe_rt0_1d
867  fe_dof_bary( 1 , 1 , f ) = 1.0_rp
868  fe_dof_bary( 1 , 2 , f ) = 1.0_rp
869 
871  x1 = 0.5_rp * ( 1.0_rp + 1.0_rp / sqrt(3.0_rp ) )
872  x2 = 0.5_rp * ( 1.0_rp - 1.0_rp / sqrt(3.0_rp ) )
873  fe_dof_bary( 1: 2, 1, f ) = (/ x1, x2 /)
874  fe_dof_bary( 1: 2, 2, f ) = (/ x2, x1 /)
875 
876  f = fe_rt1_1d
877  fe_dof_bary( 1 , 1, f ) = 1.0_rp
878  fe_dof_bary( 1 , 2, f ) = 1.0_rp
879  fe_dof_bary( 1: 2, 3, f ) = (/0.5_rp, 0.5_rp/)
880 
881  f = fe_rt1_1d_dual
882  fe_dof_bary( 1 , 1, f ) = 1.0_rp
883  fe_dof_bary( 1 , 2, f ) = 1.0_rp
884  fe_dof_bary( 1: 2, 3, f ) = (/0.5_rp, 0.5_rp/)
885 
886  f = fe_rt1_1d_2
887  fe_dof_bary( 1 , 1, f ) = 1.0_rp
888  fe_dof_bary( 1 , 2, f ) = 1.0_rp
889  fe_dof_bary( 1: 2, 3, f ) = (/0.5_rp, 0.5_rp/)
890 
891  f = fe_rt1_1d_2_dual
892  fe_dof_bary( 1 , 1, f ) = 1.0_rp
893  fe_dof_bary( 1 , 2, f ) = 1.0_rp
894  fe_dof_bary( 1: 2, 3, f ) = (/0.5_rp, 0.5_rp/)
895 
896  f = fe_rt0_1d_dual
897  fe_dof_bary( 1 , 1 , f ) = 1.0_rp
898  fe_dof_bary( 1 , 2 , f ) = 1.0_rp
899 
900  end subroutine def_fe_dof_bary
901 
904  subroutine def_fe_func()
906  !! scalar finite element functions
907  !!
908  fe_func(fe_p1_1d)%u => p1_1d_u
909  fe_func(fe_p2_1d)%u => p2_1d_u
910  fe_func(fe_p3_1d)%u => p3_1d_u
911  fe_func(fe_p0_2d)%u => p0_2d_u
912  fe_func(fe_p1_2d)%u => p1_2d_u
913  fe_func(fe_p2_2d)%u => p2_2d_u
914  fe_func(fe_p3_2d)%u => p3_2d_u
915  fe_func(fe_p1_3d)%u => p1_3d_u
916  fe_func(fe_p2_3d)%u => p2_3d_u
918  fe_func(fe_p0_1d)%u => p0_1d_u
920 
921  !! gradient of scalar finite element functions
922  !!
923  fe_func(fe_p1_1d)%grad_u => p1_1d_grad_u
924  fe_func(fe_p2_1d)%grad_u => p2_1d_grad_u
925  fe_func(fe_p3_1d)%grad_u => p3_1d_grad_u
926  fe_func(fe_p1_2d)%grad_u => p1_2d_grad_u
927  fe_func(fe_p2_2d)%grad_u => p2_2d_grad_u
928  fe_func(fe_p3_2d)%grad_u => p3_2d_grad_u
929  fe_func(fe_p1_3d)%grad_u => p1_3d_grad_u
930  fe_func(fe_p2_3d)%grad_u => p2_3d_grad_u
931 
932  !! vector finite element functions
933  !!
934  fe_func(fe_rt0_2d )%phi => rt0_2d_phi
936  fe_func(fe_rt0_1d )%phi => rt0_1d_phi
937  fe_func(fe_rt1_1d )%phi => rt1_1d_phi
942 
943  !! divergence of vector finite element functions
944  !!
945  fe_func(fe_rt0_2d)%div_phi => rt0_2d_div_phi
947  fe_func(fe_rt0_1d)%div_phi => rt0_1d_div_phi
948  fe_func(fe_rt1_1d)%div_phi => rt1_1d_div_phi
953 
954  end subroutine def_fe_func
955 
956 
959  subroutine def_fe_dof_nbvtx()
961  integer :: f, i, geo
962 
963  do f=1, fe_tot_nb
964 
965  do i=1, fe_nbdof(f)
966  geo = fe_dof_geo(i,f)
967 
968  fe_dof_nbvtx(i,f) = cell_nbvtx(geo)
969 
970  end do
971  end do
972 
973  end subroutine def_fe_dof_nbvtx
974 
975 
978  subroutine def_fe_dim()
980  integer :: ft, dim, geo
981 
982  do ft=1, fe_tot_nb
983  geo = fe_geo(ft)
984  dim = cell_dim(geo)
985  fe_dim(ft) = dim
986  end do
987 
988  end subroutine def_fe_dim
989 
990 
993  subroutine def_fe_dof_coord()
995  integer :: ft, nbDof, dim, dof_nbVtx
996  integer :: ii, dof, geo, vtx
997  real(RP), dimension(3) :: x
998 
999  do ft=1, fe_tot_nb
1000 
1001  nbdof = fe_nbdof(ft)
1002  geo = fe_geo(ft)
1003  dim = fe_dim(ft)
1004 
1005  fe_dof_coord(ft) = r_2d( dim, nbdof)
1006 
1007  do dof=1, nbdof
1008 
1009  dof_nbvtx = fe_dof_nbvtx(dof, ft)
1010 
1011  x(1:dim) = 0.0_rp
1012  do ii=1, dof_nbvtx
1013 
1014  vtx = fe_dof_vtx(ii, dof, ft)
1015 
1016  x(1:dim) = x(1:dim) + fe_dof_bary(ii, dof, ft) * &
1017  & cell_coord(geo)%y(:,vtx)
1018 
1019  end do
1020  fe_dof_coord(ft)%y(:, dof) = x(1:dim)
1021 
1022  end do
1023  end do
1024 
1025  end subroutine def_fe_dof_coord
1026 
1027 end module fe_mod
subroutine rt0_1d_dual_div_phi(val, n, x, m)
RT0 Dual ref choral/applications/PG_star_1D_poisson/sage/PG2_1D_P1_1.sage.
integer, parameter fe_p1_1d_disc_ortho
integer, parameter fe_rt0_2d
subroutine, public fe_init(b)
initialise FE_XXX arrays
Definition: fe_mod.F90:252
integer, dimension(fe_tot_nb), public fe_dim
Reference cell dimension for each fe method.
Definition: fe_mod.F90:74
type(r_2d), dimension(fe_tot_nb), public fe_dof_coord
DOF node coordinates.
Definition: fe_mod.F90:125
integer, dimension(fe_max_nbdof, fe_tot_nb), public fe_dof_geo
FE_DOF_GEO(i, j) = Geometric type of DOF i of the FE method j.
Definition: fe_mod.F90:78
subroutine rt0_1d_phi(val, n, x, m)
RT0 1D.
integer, parameter cell_tet
Tetrahedron.
BASIC TOOLS
Definition: basic_tools.F90:8
DEFINITION OF FINITE ELEMENT BASIS FUNCTIONS
Definition: fe_func_mod.f90:32
integer, dimension(cell_tot_nb), public cell_nbvtx
Number of vertexes for each cell type.
Definition: cell_mod.f90:118
subroutine p1_3d_u(val, n, x, m)
P1_3d.
integer, parameter fe_p2_2d
subroutine p1_2d_disc_ortho_u(val, n, x, m)
P1_2d_disc_ortho.
subroutine rt1_2d_2_div_phi(val, n, xx, m)
RT1_2D.
subroutine rt1_1d_dual_phi(val, n, x, m)
RT1_1D_DUAL ref = applications/PG_star_1D_poisson/sage/PG2_1D_P2_1.sage.
subroutine rt1_1d_div_phi(val, n, x, m)
RT1_1D, considered with the usual basis for the DOFs: f(0), f(1), f(1/2) ref = applications/PG_star_1...
subroutine rt0_2d_div_phi(val, n, x, m)
RT0 2D.
integer, parameter fe_dof_nrm_trace
Normal component to a face (vector finite element)
subroutine p2_3d_u(val, n, x, m)
P2_3d.
subroutine p1_1d_u(val, n, x, m)
P1_1d.
Definition: fe_func_mod.f90:87
DEFINITION OF FINITE ELEMENT METHODS
Definition: fe_mod.F90:19
integer, parameter fe_rt1_1d_dual
subroutine p1_1d_disc_ortho_u(val, n, x, m)
P1_1d_disc_ortho.
subroutine p3_2d_grad_u(val, n, x, m)
P3_2d.
subroutine def_fe_dof_geo()
DOF geometrical types.
Definition: fe_mod.F90:445
integer, parameter fe_dof_comp1
Nodal first component (vector finite element)
subroutine rt0_1d_dual_phi(val, n, x, m)
RT0 Dual ref choral/applications/PG_star_1D_poisson/sage/PG2_1D_P1_1.sage.
subroutine, public quit(message)
Stop program execution, display an error messahe.
integer, parameter fe_p1_2d_disc_ortho
subroutine rt1_1d_2_dual_div_phi(val, n, x, m)
RT1_1D_2_DUAL ref = applications/PG_star_1D_poisson/sage/PG2_1D_P2_3.sage.
integer, parameter fe_dof_flx
Interface flux (vector finite element)
integer, parameter fe_rt0_1d_dual
subroutine rt1_2d_2_phi(val, n, xx, m)
RT1_2D.
integer, parameter fe_dof_comp2
Nodal second component (vector finite element)
integer, parameter fe_rt1_2d_2
integer, dimension(fe_max_nbdof, fe_tot_nb), public fe_dof_nbvtx
FE_DOF_NBVTX(i,j) = number of vertexes of the geometric element associated with the DOF i of the FE m...
Definition: fe_mod.F90:104
integer, parameter cell_edg
Edge (line segment)
integer, parameter fe_tot_nb
Number of FE methods.
subroutine p1_1d_grad_u(val, n, x, m)
P1_1d.
integer, dimension(cell_tot_nb), public cell_dim
Dimension for each cell type.
Definition: cell_mod.f90:130
conversion integers or rational to real
Definition: real_type.F90:153
integer, parameter fe_max_nbdof
Maximum number of dof for an element.
REAL NUMBERS PRECISION IN CHORAL: selects simple/double/quad
Definition: real_type.F90:33
subroutine, public scalgrad_fe(val, nbDof, x, dim, ft)
Evaluate the gradient of the FE basis functions at point x.
Definition: fe_mod.F90:166
integer, parameter fe_p2_3d
subroutine p2_2d_u(val, n, x, m)
P2_2d.
subroutine p0_2d_u(val, n, x, m)
P0_2d.
subroutine def_fe()
Number of DOF for each FE method.
Definition: fe_mod.F90:276
subroutine def_fe_dim()
set FE dimension
Definition: fe_mod.F90:979
subroutine p2_3d_grad_u(val, n, x, m)
P2_3d.
subroutine rt1_1d_phi(val, n, x, m)
RT1_1D, considered with the usual basis for the DOFs: f(0), f(1), f(1/2) ref = applications/PG_star_1...
integer, parameter fe_p3_1d
type(r_2d), dimension(cell_tot_nb), public cell_coord
Cell node coordinates.
Definition: cell_mod.f90:155
subroutine p2_1d_grad_u(val, n, x, m)
P2_1d.
integer, parameter cell_trg
Triangle.
integer, dimension(fe_tot_nb), public fe_nbdof
Number of DOF for each fe method.
Definition: fe_mod.F90:68
subroutine p1_2d_grad_u(val, n, x, m)
P1_2d.
subroutine def_fe_dof_geoelmt()
in the cell K :
Definition: fe_mod.F90:549
CHORAL CONSTANTS
subroutine, public vect_fe(val, nbDof, x, dim, ft)
Evaluate vector FE basis functions at point x.
Definition: fe_mod.F90:195
subroutine p0_1d_u(val, n, x, m)
P0_1d.
subroutine def_fe_dof_vtx()
FE_DOF_VTX(1:n, i, f) = vertexes of the ref. cell associated to the dof i of the FE method f...
Definition: fe_mod.F90:642
subroutine p3_2d_u(val, n, x, m)
P3_2d.
type(fe_functions), dimension(fe_tot_nb), public fe_func
FE_ARRAY Arrays describing finite element methods.
Definition: fe_mod.F90:62
subroutine rt1_1d_2_phi(val, n, x, m)
RT1_1D_2, considered with the alternative basis for the DOFs: f(0), f(1), ^1 f(x) dx ref = applicatio...
subroutine p1_2d_u(val, n, x, m)
P1_2d.
integer, parameter fe_p1_2d
subroutine rt0_2d_phi(val, n, x, m)
RT0 2D.
subroutine rt1_1d_2_div_phi(val, n, x, m)
RT1 1D_2, considered with the alternative basis for the DOFs: f(0), f(1), ^1 f(x) dx ref = applicatio...
subroutine p2_2d_grad_u(val, n, x, m)
P2_2d.
integer, parameter fe_rt1_1d_2_dual
integer, dimension(cell_max_nbvtx, fe_max_nbdof, fe_tot_nb), public fe_dof_vtx
FE_DOF_VTX(1:n,i,j) = for the DOF i of the FE method j :
Definition: fe_mod.F90:112
subroutine rt1_1d_2_dual_phi(val, n, x, m)
RT1_1D_2_DUAL ref = applications/PG_star_1D_poisson/sage/PG2_1D_P2_3.sage.
subroutine p1_3d_grad_u(val, n, x, m)
P1_3d.
integer, parameter fe_dof_lag
Nodal value (Lagrangian DOF, scalar finite element)
integer, dimension(fe_tot_nb), public fe_geo
Reference cell geometry for each fe method.
Definition: fe_mod.F90:71
character(len=16), dimension(0:fe_tot_nb), public fe_name
Name for each fe method.
Definition: fe_mod.F90:65
subroutine def_fe_dof_bary()
FE_DOF_BARY(1:n, i, f) = barycentric coordinates of the dof i of the FE method f. ...
Definition: fe_mod.F90:776
integer, parameter fe_rt1_1d
subroutine rt1_1d_dual_div_phi(val, n, x, m)
RT1_1D_DUAL ref = applications/PG_star_1D_poisson/sage/PG2_1D_P2_1.sage.
subroutine def_fe_dof_type()
DOF types.
Definition: fe_mod.F90:350
integer, parameter fe_none
subroutine rt0_1d_div_phi(val, n, x, m)
RT0 1D.
integer, parameter fe_rt1_1d_2
subroutine p2_1d_u(val, n, x, m)
P2_1d.
Definition: fe_func_mod.f90:98
integer, parameter fe_p0_1d
integer, parameter cell_vtx
Vertex.
integer, dimension(fe_max_nbdof, fe_tot_nb), public fe_dof_geoelmt
For the FE method j : the associated reference cell is := FE_GEO(i)
Definition: fe_mod.F90:99
subroutine def_fe_dof_coord()
set FE dof coordinates
Definition: fe_mod.F90:994
real(rp), dimension(cell_max_nbvtx, fe_max_nbdof, fe_tot_nb), public fe_dof_bary
FE_DOF_BARY(1:n,i,j) = for the DOF i of the FE method j :
Definition: fe_mod.F90:121
integer, parameter fe_p3_2d
subroutine, public vectdiv_fe(val, nbDof, x, dim, ft)
Evaluate the divergence of the vector FE basis functions at point x.
Definition: fe_mod.F90:225
R_2D: this type will allow to define arrays of 2D real arrays.
Definition: real_type.F90:114
integer, parameter fe_p0_2d
integer, dimension(fe_max_nbdof, fe_tot_nb), public fe_dof_type
FE_DOF_TYPE(i,j) = dof type of the DOF i of the FE method j (Lagrange, ...)
Definition: fe_mod.F90:82
integer, parameter fe_p1_3d
integer, parameter fe_p1_1d
Functions associated with a finite element.
Definition: fe_func_mod.f90:41
subroutine def_fe_func()
Finite element function.
Definition: fe_mod.F90:905
DEFINITION OF GEOMETRICAL CELLS (for meshes)
Definition: cell_mod.f90:82
subroutine p3_1d_u(val, n, x, m)
P3_1d.
subroutine def_fe_dof_nbvtx()
FE DOF number of vertexes.
Definition: fe_mod.F90:960
integer, parameter fe_rt0_1d
integer, parameter fe_p2_1d
subroutine, public scal_fe(val, nbDof, x, dim, ft)
Evaluate FE basis functions at point x.
Definition: fe_mod.F90:136
subroutine p3_1d_grad_u(val, n, x, m)
P3_1d.