Choral
feSpace_visu.f90
Go to the documentation of this file.
1 
11 
12 program fespace_visu
13 
14  use choral
16 
17  implicit none
18 
19  !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
20  !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
21  !!
22  !! VARAIBLE DEFINITION : BEGIN
23  !!
24  !!
25 
26  !! verb = verbosity level
27  !!
28  integer, parameter :: verb = 2
29 
30  !! fe = finite element type
31  !!
32  integer, parameter :: fe = fe_p2_2d
33 
34  !! msh_file = mesh file name
35  !! msh = mesh
36  !! X_h = finite element space
37  !!
38  character(len=150) :: msh_file
39  type(mesh) :: msh
40  type(fespace) :: X_h
41 
42 
43  !! Variable to store a finite element function
44  !!
45  real(RP), dimension(:), allocatable :: f_h
46 
47  !! output files
48  !!
49  character(len=150) :: feSpace_file, str
50  !!
51  !! VARAIBLE DEFINITION : END
52  !!
53  !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
54  !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
55 
56  call choral_init(verb=verb)
57 
58  write(*,*)
59  write(*,*)'mesh_visu: start'
60  write(*,*)""
61  write(*,*)" ASSEMBLE A MESH "
62  write(*,*)""
63  write(*,*)" and visualise it with gmsh"
64  write(*,*)""
65 
66  !!
67  !! settings: choose a mesh
68  !!
69  ! msh_file = trim(GMSH_DIR)//"disk/disk1_1.msh"
70  msh_file = trim(gmsh_dir)//"testMesh/square.msh"
71  ! msh_file = trim(GMSH_DIR)//"sphere/sphere1_1.msh"
72  ! msh_file = trim(GMSH_DIR)//"sphere/sphere2_1.msh"
73 
74 
75  write(*,*) ""
76  write(*,*)"================================== ASSEMBLING THE MESH "
77  msh = mesh(trim(msh_file), "gmsh")
78  if (verb>1) call print(msh)
79 
80 
81  write(*,*) ""
82  write(*,*)"================================== ASSEMBLING THE&
83  & FINITE ELEMENT SPACE X_h"
84  x_h = fespace(msh)
85  call set(x_h, fe)
86  call assemble(x_h)
87 
88 
89  write(*,*) ""
90  write(*,*)"================================== DISPLAY THE&
91  & FINITE ELEMENT NODES"
92  fespace_file = "feSpace.msh"
93  call write(x_h, fespace_file, "gmsh")
94  call system("gmsh -option &
95  &"//trim(gmsh_dir)//"gmsh-options-mesh "//trim(fespace_file))
96 
97 
98  write(*,*) ""
99  write(*,*)"================================== DISPLAY A&
100  & FUNCTION INTERPOLATED ON X_h"
101  !!
102  !! 1. Interpolate 'f' to the finite element function f_h \in X_h
103  call interp_scal_func(f_h, f, x_h)
104  !!
105  !! 2. Add the finite element function 'f_h' to the file 'feSpace_file'
106  call gmsh_addview(x_h, f_h, fespace_file, &
107  & 'Interpolant of the Function f', 0._rp, 0)
108  !!
109  !! 3. Visualise
110  str='gmsh -option '//trim(gmsh_dir)//'gmsh-options-view '&
111  &//trim(fespace_file)
112  call system(trim(str))
113 
114 
115 contains
116 
117  function f(x)
118  real(RP), dimension(3), intent(in) :: x
119  real(RP) :: f
120 
121  f = sin(2*pi*x(1)) * sin(2*pi*x(2))
122 
123  end function f
124 
125 
126 
127 end program fespace_visu
subroutine, public choral_init(verb)
Initialisation for the CHORAL library
Definition: choral.F90:166
program fespace_visu
ASSEMBLE A FINITE ELEMENS SPACE of type feSpace
TOP-LEVEL MODULE FOR THE LIBRARY CHORAL
Definition: choral.F90:27
integer, parameter fe_p2_2d
CHORAL CONSTANTS
real(rp) function f(x)
right hand side 'f'