xref: /libCEED/examples/ceed/ex1-volume-f.f90 (revision d6ed5abffa3823f6603d4dc65272439584c38270)
1*860dc821SJeremy L Thompson! Copyright (c) 2017-2026,  Lawrence Livermore National Security,  LLC and other CEED contributors.
2*860dc821SJeremy L Thompson! All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*860dc821SJeremy L Thompson
4*860dc821SJeremy L Thompson! SPDX-License-Identifier: BSD-2-Clause
5*860dc821SJeremy L Thompson
6*860dc821SJeremy L Thompson! This file is part of CEED:  http:Cgithub.com/ceed
7*860dc821SJeremy L Thompson
8*860dc821SJeremy L Thompson! libCEED Example 1
9*860dc821SJeremy L Thompson
10*860dc821SJeremy L Thompson! This example illustrates a simple usage of libCEED to compute the volume of a 3D body using matrix-free application of a mass operator.
11*860dc821SJeremy L Thompson! Arbitrary mesh and solution degrees in 1D,  2D and 3D are supported from the same code.
12*860dc821SJeremy L Thompson
13*860dc821SJeremy L Thompson! The example has no dependencies,  and is designed to be self-contained.
14*860dc821SJeremy L Thompson! For additional examples that use external discretization libraries (MFEM,  PETSc,  etc.) see the subdirectories in libceed/examples.
15*860dc821SJeremy L Thompson
16*860dc821SJeremy L Thompson! All libCEED objects use a Ceed device object constructed based on a command line argument (-ceed).
17*860dc821SJeremy L Thompson
18*860dc821SJeremy L Thompson! Build with:
19*860dc821SJeremy L Thompson
20*860dc821SJeremy L Thompson!     make ex1-volume [CEED_DIR = </path/to/libceed>]
21*860dc821SJeremy L Thompson
22*860dc821SJeremy L Thompson! Sample runs:
23*860dc821SJeremy L Thompson
24*860dc821SJeremy L Thompson!     ./ex1-volume-f
25*860dc821SJeremy L Thompson!     ./ex1-volume-f -ceed /cpu/self
26*860dc821SJeremy L Thompson!     ./ex1-volume-f -ceed /gpu/cuda
27*860dc821SJeremy L Thompson
28*860dc821SJeremy L Thompson! Test in 1D-3D
29*860dc821SJeremy L Thompson! TESTARGS(name = "1D User QFunction") -ceed {ceed_resource} -d 1 -t
30*860dc821SJeremy L Thompson! TESTARGS(name = "2D User QFunction") -ceed {ceed_resource} -d 2 -t
31*860dc821SJeremy L Thompson! TESTARGS(name = "3D User QFunction") -ceed {ceed_resource} -d 3 -t
32*860dc821SJeremy L Thompson! TESTARGS(name = "1D Gallery QFunction") -ceed {ceed_resource} -d 1 -t -g
33*860dc821SJeremy L Thompson! TESTARGS(name = "2D Gallery QFunction") -ceed {ceed_resource} -d 2 -t -g
34*860dc821SJeremy L Thompson! TESTARGS(name = "3D Gallery QFunction") -ceed {ceed_resource} -d 3 -t -g
35*860dc821SJeremy L Thompson
36*860dc821SJeremy L Thompson!> @file
37*860dc821SJeremy L Thompson!> libCEED example using mass operator to compute volume
38*860dc821SJeremy L Thompson
39*860dc821SJeremy L Thompson    include 'ex1-volume-f.h'
40*860dc821SJeremy L Thompson
41*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
42*860dc821SJeremy L Thompsonsubroutine getcartesianmeshsize(fe_dim, degree, prob_size, num_xyz)
43*860dc821SJeremy L Thompson    implicit none
44*860dc821SJeremy L Thompson    integer fe_dim
45*860dc821SJeremy L Thompson    integer degree
46*860dc821SJeremy L Thompson    integer prob_size
47*860dc821SJeremy L Thompson    integer num_xyz(3)
48*860dc821SJeremy L Thompson
49*860dc821SJeremy L Thompson    integer num_elem
50*860dc821SJeremy L Thompson    integer s, r, d, sd
51*860dc821SJeremy L Thompson    num_elem = prob_size/(degree**fe_dim)
52*860dc821SJeremy L Thompson    s = 0
53*860dc821SJeremy L Thompson
54*860dc821SJeremy L Thompson! Use the approximate formula:
55*860dc821SJeremy L Thompson!    prob_size ~ num_elem * degree^dim
56*860dc821SJeremy L Thompson! find s: num_elem/2 < 2^s <= num_elem
57*860dc821SJeremy L Thompson
58*860dc821SJeremy L Thompson  do while (num_elem > 1)
59*860dc821SJeremy L Thompson    num_elem = num_elem/2
60*860dc821SJeremy L Thompson    s = s + 1
61*860dc821SJeremy L Thompson  end do
62*860dc821SJeremy L Thompson  r = mod(s, fe_dim)
63*860dc821SJeremy L Thompson
64*860dc821SJeremy L Thompson  do d = 1, fe_dim
65*860dc821SJeremy L Thompson    sd = s/fe_dim
66*860dc821SJeremy L Thompson    if (r > 0) then
67*860dc821SJeremy L Thompson      sd = sd + 1
68*860dc821SJeremy L Thompson      r = r - 1
69*860dc821SJeremy L Thompson    end if
70*860dc821SJeremy L Thompson    num_xyz(d) = ISHFT(1, sd)
71*860dc821SJeremy L Thompson  end do
72*860dc821SJeremy L Thompsonend
73*860dc821SJeremy L Thompson
74*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
75*860dc821SJeremy L Thompsonsubroutine buildcartesianrestriction(ceed, fe_dim, num_xyz, degree, num_comp, mesh_size, num_qpts, restriction,&
76*860dc821SJeremy L Thompson&     q_data_restriction, err)
77*860dc821SJeremy L Thompson    implicit none
78*860dc821SJeremy L Thompson    include 'ceed/fortran.h'
79*860dc821SJeremy L Thompson
80*860dc821SJeremy L Thompson    integer ceed
81*860dc821SJeremy L Thompson    integer fe_dim
82*860dc821SJeremy L Thompson    integer num_xyz(3)
83*860dc821SJeremy L Thompson    integer degree
84*860dc821SJeremy L Thompson    integer num_comp
85*860dc821SJeremy L Thompson    integer mesh_size
86*860dc821SJeremy L Thompson    integer num_qpts
87*860dc821SJeremy L Thompson    integer restriction
88*860dc821SJeremy L Thompson    integer q_data_restriction
89*860dc821SJeremy L Thompson    integer err
90*860dc821SJeremy L Thompson
91*860dc821SJeremy L Thompson    integer p
92*860dc821SJeremy L Thompson    integer num_nodes
93*860dc821SJeremy L Thompson    integer elem_qpts
94*860dc821SJeremy L Thompson    integer num_elem
95*860dc821SJeremy L Thompson    integer scalar_size
96*860dc821SJeremy L Thompson    integer nd(3)
97*860dc821SJeremy L Thompson    integer elem_nodes_size
98*860dc821SJeremy L Thompson    integer e_xyz(3),  re
99*860dc821SJeremy L Thompson    integer g_nodes, g_nodes_stride, r_nodes
100*860dc821SJeremy L Thompson    integer, dimension (:), allocatable :: elem_nodes
101*860dc821SJeremy L Thompson
102*860dc821SJeremy L Thompson    integer i, j, k
103*860dc821SJeremy L Thompson
104*860dc821SJeremy L Thompson    p = degree + 1
105*860dc821SJeremy L Thompson    num_nodes = p**fe_dim
106*860dc821SJeremy L Thompson    elem_qpts = num_qpts**fe_dim
107*860dc821SJeremy L Thompson    num_elem  = 1
108*860dc821SJeremy L Thompson    scalar_size = 1
109*860dc821SJeremy L Thompson
110*860dc821SJeremy L Thompson    do i = 1, fe_dim
111*860dc821SJeremy L Thompson      num_elem = num_elem * num_xyz(i)
112*860dc821SJeremy L Thompson      nd(i) = num_xyz(i) * (p - 1) + 1
113*860dc821SJeremy L Thompson      scalar_size = scalar_size*nd(i)
114*860dc821SJeremy L Thompson    end do
115*860dc821SJeremy L Thompson    mesh_size = scalar_size*num_comp
116*860dc821SJeremy L Thompson! elem:       0         1             n-1
117*860dc821SJeremy L Thompson!         |---*-...-*---|---*-...-*---|- ... -|--...--|
118*860dc821SJeremy L Thompson! num_nodes:   0   1    p-1  p  p+1     2*p         n*p
119*860dc821SJeremy L Thompson    elem_nodes_size = num_elem*num_nodes
120*860dc821SJeremy L Thompson    allocate (elem_nodes(elem_nodes_size))
121*860dc821SJeremy L Thompson
122*860dc821SJeremy L Thompson    do i = 1, num_elem
123*860dc821SJeremy L Thompson      e_xyz(1) = 1
124*860dc821SJeremy L Thompson      e_xyz(2) = 1
125*860dc821SJeremy L Thompson      e_xyz(3) = 1
126*860dc821SJeremy L Thompson      re = i - 1
127*860dc821SJeremy L Thompson
128*860dc821SJeremy L Thompson      do j = 1, fe_dim
129*860dc821SJeremy L Thompson        e_xyz(j) = mod(re, num_xyz(j))
130*860dc821SJeremy L Thompson        re = re/num_xyz(j)
131*860dc821SJeremy L Thompson      end do
132*860dc821SJeremy L Thompson
133*860dc821SJeremy L Thompson      do j = 1, num_nodes
134*860dc821SJeremy L Thompson        g_nodes = 0
135*860dc821SJeremy L Thompson        g_nodes_stride = 1
136*860dc821SJeremy L Thompson        r_nodes = j - 1
137*860dc821SJeremy L Thompson
138*860dc821SJeremy L Thompson        do k = 1, fe_dim
139*860dc821SJeremy L Thompson          g_nodes = g_nodes + (e_xyz(k) * (p - 1) + mod(r_nodes, p)) * g_nodes_stride
140*860dc821SJeremy L Thompson          g_nodes_stride = g_nodes_stride * nd(k)
141*860dc821SJeremy L Thompson          r_nodes = r_nodes/p
142*860dc821SJeremy L Thompson        end do
143*860dc821SJeremy L Thompson        elem_nodes((i - 1) * num_nodes + j) = g_nodes
144*860dc821SJeremy L Thompson      end do
145*860dc821SJeremy L Thompson    end do
146*860dc821SJeremy L Thompson
147*860dc821SJeremy L Thompson    call ceedelemrestrictioncreate(ceed, num_elem, num_nodes, num_comp, scalar_size, mesh_size, ceed_mem_host,&
148*860dc821SJeremy L Thompson             &ceed_copy_values, elem_nodes, restriction, err)
149*860dc821SJeremy L Thompson    if (q_data_restriction /=  ceed_qfunction_none) then
150*860dc821SJeremy L Thompson      call ceedelemrestrictioncreatestrided(ceed, num_elem, elem_qpts, num_comp, num_comp * elem_qpts * num_elem,&
151*860dc821SJeremy L Thompson               &ceed_strides_backend, q_data_restriction, err)
152*860dc821SJeremy L Thompson    end if
153*860dc821SJeremy L Thompson    deallocate (elem_nodes)
154*860dc821SJeremy L Thompsonend
155*860dc821SJeremy L Thompson
156*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
157*860dc821SJeremy L Thompsonsubroutine transformmeshcoords(fe_dim, mesh_size, coords, exact_volume, err)
158*860dc821SJeremy L Thompson    implicit none
159*860dc821SJeremy L Thompson
160*860dc821SJeremy L Thompson    integer fe_dim
161*860dc821SJeremy L Thompson    integer mesh_size, scalar_size
162*860dc821SJeremy L Thompson    real*8 coords(mesh_size)
163*860dc821SJeremy L Thompson    real*8 exact_volume
164*860dc821SJeremy L Thompson    real*8 m_pi, m_pi_2
165*860dc821SJeremy L Thompson    parameter(m_pi = 3.14159265358979323846d0)
166*860dc821SJeremy L Thompson    parameter(m_pi_2 = 1.57079632679489661923d0)
167*860dc821SJeremy L Thompson    integer err
168*860dc821SJeremy L Thompson
169*860dc821SJeremy L Thompson    integer i
170*860dc821SJeremy L Thompson    real*8 u, v
171*860dc821SJeremy L Thompson
172*860dc821SJeremy L Thompson    scalar_size = mesh_size/fe_dim
173*860dc821SJeremy L Thompson    select case (fe_dim)
174*860dc821SJeremy L Thompson    case (1)
175*860dc821SJeremy L Thompson      do i = 1, scalar_size
176*860dc821SJeremy L Thompson        coords(i) = 0.5d0 + (1.d0/sqrt(3.d0)) * sin((2.d0/3.d0) * m_pi * (coords(i) - 0.5d0))
177*860dc821SJeremy L Thompson      end do
178*860dc821SJeremy L Thompson      exact_volume = 1.d0
179*860dc821SJeremy L Thompson
180*860dc821SJeremy L Thompson    case (2,  3)
181*860dc821SJeremy L Thompson      do i = 1, scalar_size
182*860dc821SJeremy L Thompson        u = 1.d0 + coords(i)
183*860dc821SJeremy L Thompson        v = m_pi_2 * coords(i + scalar_size)
184*860dc821SJeremy L Thompson
185*860dc821SJeremy L Thompson        coords(i)               = u * cos(v)
186*860dc821SJeremy L Thompson        coords(i + scalar_size) = u * sin(v)
187*860dc821SJeremy L Thompson      end do
188*860dc821SJeremy L Thompson      exact_volume = 3.d0/4.d0 * m_pi
189*860dc821SJeremy L Thompson    end select
190*860dc821SJeremy L Thompsonend
191*860dc821SJeremy L Thompson
192*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
193*860dc821SJeremy L Thompsonsubroutine setcartesianmeshcoords(fe_dim, num_xyz, mesh_degree, mesh_coords, exact_volume, err)
194*860dc821SJeremy L Thompson    implicit none
195*860dc821SJeremy L Thompson    include 'ceed/fortran.h'
196*860dc821SJeremy L Thompson
197*860dc821SJeremy L Thompson    integer fe_dim
198*860dc821SJeremy L Thompson    integer num_xyz(3)
199*860dc821SJeremy L Thompson    integer mesh_degree
200*860dc821SJeremy L Thompson    integer mesh_coords
201*860dc821SJeremy L Thompson    real*8 exact_volume
202*860dc821SJeremy L Thompson    integer err
203*860dc821SJeremy L Thompson
204*860dc821SJeremy L Thompson    integer p
205*860dc821SJeremy L Thompson    integer scalar_size
206*860dc821SJeremy L Thompson    integer coords_size
207*860dc821SJeremy L Thompson    integer r_nodes
208*860dc821SJeremy L Thompson    integer d_1d
209*860dc821SJeremy L Thompson    integer nd(3)
210*860dc821SJeremy L Thompson    real*8, dimension (:), allocatable :: nodes,  qpts
211*860dc821SJeremy L Thompson    real*8, dimension (:), allocatable :: coords
212*860dc821SJeremy L Thompson    integer*8 offset
213*860dc821SJeremy L Thompson    integer i, j
214*860dc821SJeremy L Thompson    p = mesh_degree + 1
215*860dc821SJeremy L Thompson    scalar_size = 1
216*860dc821SJeremy L Thompson
217*860dc821SJeremy L Thompson    do i = 1, fe_dim
218*860dc821SJeremy L Thompson      nd(i) = num_xyz(i) * (p - 1) + 1
219*860dc821SJeremy L Thompson      scalar_size = scalar_size * nd(i)
220*860dc821SJeremy L Thompson    end do
221*860dc821SJeremy L Thompson
222*860dc821SJeremy L Thompson    coords_size = scalar_size * fe_dim
223*860dc821SJeremy L Thompson    allocate (coords(coords_size))
224*860dc821SJeremy L Thompson
225*860dc821SJeremy L Thompson! The H1 basis uses Lobatto quadrature points as nodes
226*860dc821SJeremy L Thompson    allocate (nodes(p))
227*860dc821SJeremy L Thompson    allocate (qpts(p))
228*860dc821SJeremy L Thompson    call ceedlobattoquadrature(p, nodes, qpts, err)
229*860dc821SJeremy L Thompson    deallocate(qpts)
230*860dc821SJeremy L Thompson    do i = 1, p
231*860dc821SJeremy L Thompson      nodes(i) = 0.5 + 0.5 * nodes(i)
232*860dc821SJeremy L Thompson    end do
233*860dc821SJeremy L Thompson
234*860dc821SJeremy L Thompson    do i = 1, scalar_size
235*860dc821SJeremy L Thompson      r_nodes = i - 1
236*860dc821SJeremy L Thompson
237*860dc821SJeremy L Thompson      do j = 1, fe_dim
238*860dc821SJeremy L Thompson        d_1d  =  mod(r_nodes, nd(j))
239*860dc821SJeremy L Thompson        coords(scalar_size * (j - 1) + i) = ((d_1d/(p - 1)) + nodes(mod(d_1d, p - 1) + 1))/num_xyz(j)
240*860dc821SJeremy L Thompson        r_nodes = r_nodes/nd(j)
241*860dc821SJeremy L Thompson      end do
242*860dc821SJeremy L Thompson    end do
243*860dc821SJeremy L Thompson    deallocate(nodes)
244*860dc821SJeremy L Thompson
245*860dc821SJeremy L Thompson    call transformmeshcoords(fe_dim, coords_size, coords, exact_volume, err)
246*860dc821SJeremy L Thompson
247*860dc821SJeremy L Thompson    offset = 0
248*860dc821SJeremy L Thompson    call ceedvectorsetarray(mesh_coords, ceed_mem_host, ceed_copy_values, coords, offset, err)
249*860dc821SJeremy L Thompson    deallocate(coords)
250*860dc821SJeremy L Thompsonend
251*860dc821SJeremy L Thompson
252*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
253*860dc821SJeremy L Thompsonprogram main
254*860dc821SJeremy L Thompson    implicit none
255*860dc821SJeremy L Thompson    include 'ceed/fortran.h'
256*860dc821SJeremy L Thompson
257*860dc821SJeremy L Thompson    character ceed_spec*32
258*860dc821SJeremy L Thompson    integer fe_dim, num_comp_x, mesh_degree, sol_degree, num_qpts
259*860dc821SJeremy L Thompson    integer num_elem, num_xyz(3), elem_qpts
260*860dc821SJeremy L Thompson    integer prob_size, mesh_size, sol_size
261*860dc821SJeremy L Thompson    integer help, test, gallery, benchmark
262*860dc821SJeremy L Thompson    integer i, num_args, err
263*860dc821SJeremy L Thompson    character arg*32, arg_value*32
264*860dc821SJeremy L Thompson    real*8 exact_volume, computed_volume
265*860dc821SJeremy L Thompson
266*860dc821SJeremy L Thompson    integer ceed
267*860dc821SJeremy L Thompson    real*8, dimension (:), allocatable :: u_array, v_array
268*860dc821SJeremy L Thompson    integer mesh_coords, q_data, u, v
269*860dc821SJeremy L Thompson    integer mesh_restriction, sol_restriction, q_data_restriction
270*860dc821SJeremy L Thompson    integer mesh_basis, sol_basis
271*860dc821SJeremy L Thompson    integer*8 offset
272*860dc821SJeremy L Thompson    integer build_ctx
273*860dc821SJeremy L Thompson    integer build_ctx_size
274*860dc821SJeremy L Thompson    parameter(build_ctx_size = 2)
275*860dc821SJeremy L Thompson    integer*8 build_ctx_data(build_ctx_size)
276*860dc821SJeremy L Thompson    integer qf_build, qf_apply
277*860dc821SJeremy L Thompson    integer op_build, op_apply
278*860dc821SJeremy L Thompson
279*860dc821SJeremy L Thompson    external build_mass, apply_mass
280*860dc821SJeremy L Thompson
281*860dc821SJeremy L Thompson! Initial values
282*860dc821SJeremy L Thompson    ceed_spec   = '/cpu/self'
283*860dc821SJeremy L Thompson    fe_dim      = 3
284*860dc821SJeremy L Thompson    num_comp_x  = 3
285*860dc821SJeremy L Thompson    mesh_degree = 4
286*860dc821SJeremy L Thompson    sol_degree  = 4
287*860dc821SJeremy L Thompson    num_qpts    = mesh_degree + 2
288*860dc821SJeremy L Thompson    prob_size   = -1
289*860dc821SJeremy L Thompson    help      = 0
290*860dc821SJeremy L Thompson    test      = 0
291*860dc821SJeremy L Thompson    gallery   = 0
292*860dc821SJeremy L Thompson    benchmark = 0
293*860dc821SJeremy L Thompson
294*860dc821SJeremy L Thompson! Process command line arguments
295*860dc821SJeremy L Thompson
296*860dc821SJeremy L Thompson    num_args = command_argument_count()
297*860dc821SJeremy L Thompson    do i = 1, num_args
298*860dc821SJeremy L Thompson      call get_command_argument(i, arg)
299*860dc821SJeremy L Thompson
300*860dc821SJeremy L Thompson      select case (arg)
301*860dc821SJeremy L Thompson! LCOV_EXCL_START
302*860dc821SJeremy L Thompson        case ('-h')
303*860dc821SJeremy L Thompson          help = 1
304*860dc821SJeremy L Thompson
305*860dc821SJeremy L Thompson        case ('-c',  '-ceed')
306*860dc821SJeremy L Thompson          call get_command_argument(i + 1, ceed_spec)
307*860dc821SJeremy L Thompson
308*860dc821SJeremy L Thompson        case ('-d')
309*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
310*860dc821SJeremy L Thompson          read(arg_value, '(I10)') fe_dim
311*860dc821SJeremy L Thompson          num_comp_x = fe_dim
312*860dc821SJeremy L Thompson
313*860dc821SJeremy L Thompson        case ('-m')
314*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
315*860dc821SJeremy L Thompson          read(arg_value, '(I10)') mesh_degree
316*860dc821SJeremy L Thompson
317*860dc821SJeremy L Thompson        case ('-p')
318*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
319*860dc821SJeremy L Thompson          read(arg_value, '(I10)') sol_degree
320*860dc821SJeremy L Thompson
321*860dc821SJeremy L Thompson        case ('-q')
322*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
323*860dc821SJeremy L Thompson          read(arg_value, '(I10)') num_qpts
324*860dc821SJeremy L Thompson
325*860dc821SJeremy L Thompson        case ('-s')
326*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
327*860dc821SJeremy L Thompson          read(arg_value, '(I10)') prob_size
328*860dc821SJeremy L Thompson
329*860dc821SJeremy L Thompson        case ('-b')
330*860dc821SJeremy L Thompson          call get_command_argument(i + 1, arg_value)
331*860dc821SJeremy L Thompson          read(arg_value, '(I10)') benchmark
332*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
333*860dc821SJeremy L Thompson
334*860dc821SJeremy L Thompson        case ('-t')
335*860dc821SJeremy L Thompson          test = 1
336*860dc821SJeremy L Thompson
337*860dc821SJeremy L Thompson        case ('-g')
338*860dc821SJeremy L Thompson          gallery = 1
339*860dc821SJeremy L Thompson      end select
340*860dc821SJeremy L Thompson    end do
341*860dc821SJeremy L Thompson
342*860dc821SJeremy L Thompson    if (prob_size < 0) then
343*860dc821SJeremy L Thompson      if (test == 1) then
344*860dc821SJeremy L Thompson        prob_size = 8 * 16
345*860dc821SJeremy L Thompson      else
346*860dc821SJeremy L Thompson        prob_size = 256 * 1024
347*860dc821SJeremy L Thompson      end if
348*860dc821SJeremy L Thompson    end if
349*860dc821SJeremy L Thompson
350*860dc821SJeremy L Thompson! Print options
351*860dc821SJeremy L Thompson    if ((test /= 1) .OR. (help == 1)) then
352*860dc821SJeremy L Thompson! LCOV_EXCL_START
353*860dc821SJeremy L Thompson    write (*, *) 'Selected options: [command line option] : <current value>'
354*860dc821SJeremy L Thompson    write (*, *) '  Ceed specification     [-c] : ', ceed_spec
355*860dc821SJeremy L Thompson    write (*, *) '  Mesh dimension         [-d] : ', fe_dim
356*860dc821SJeremy L Thompson    write (*, *) '  Mesh degree            [-m] : ', mesh_degree
357*860dc821SJeremy L Thompson    write (*, *) '  Solution degree        [-p] : ', sol_degree
358*860dc821SJeremy L Thompson    write (*, *) '  Num. 1D quadrature pts [-q] : ', num_qpts
359*860dc821SJeremy L Thompson    write (*, *) '  Approx. # unknowns     [-s] : ', prob_size
360*860dc821SJeremy L Thompson    if (gallery == 1) then
361*860dc821SJeremy L Thompson      write (*, *) '  QFunction source       [-g] : gallery'
362*860dc821SJeremy L Thompson    else
363*860dc821SJeremy L Thompson      write (*, *) '  QFunction source       [-g] : header'
364*860dc821SJeremy L Thompson    end if
365*860dc821SJeremy L Thompson    if (help == 1) then
366*860dc821SJeremy L Thompson      if (test == 0) then
367*860dc821SJeremy L Thompson        write (*, *) 'Test/quiet mode is OFF (use -t to enable)'
368*860dc821SJeremy L Thompson      else
369*860dc821SJeremy L Thompson        write (*, *) 'Test/quiet mode is ON'
370*860dc821SJeremy L Thompson      end if
371*860dc821SJeremy L Thompson    end if
372*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
373*860dc821SJeremy L Thompson    end if
374*860dc821SJeremy L Thompson
375*860dc821SJeremy L Thompson! Select appropriate backend and logical device based on the (-ceed) command line argument
376*860dc821SJeremy L Thompson    call ceedinit(trim(ceed_spec)//char(0), ceed, err)
377*860dc821SJeremy L Thompson
378*860dc821SJeremy L Thompson! Construct the mesh and solution bases
379*860dc821SJeremy L Thompson    call ceedbasiscreatetensorh1lagrange(ceed, fe_dim, num_comp_x, mesh_degree + 1, num_qpts, ceed_gauss, mesh_basis,&
380*860dc821SJeremy L Thompson             &err)
381*860dc821SJeremy L Thompson    call ceedbasiscreatetensorh1lagrange(ceed, fe_dim, 1, sol_degree + 1, num_qpts, ceed_gauss, sol_basis, err)
382*860dc821SJeremy L Thompson
383*860dc821SJeremy L Thompson! Determine the mesh size based on the given approximate problem size
384*860dc821SJeremy L Thompson    call getcartesianmeshsize(fe_dim, sol_degree, prob_size, num_xyz)
385*860dc821SJeremy L Thompson    if (test == 0) then
386*860dc821SJeremy L Thompson! LCOV_EXCL_START
387*860dc821SJeremy L Thompson    write (*, '(A16, I8)', advance='no') 'Mesh size: nx = ', num_xyz(1)
388*860dc821SJeremy L Thompson    if (num_comp_x > 1) then
389*860dc821SJeremy L Thompson      write (*, '(A7, I8)', advance='no') ',  ny = ', num_xyz(2)
390*860dc821SJeremy L Thompson    end if
391*860dc821SJeremy L Thompson    if (num_comp_x > 2) then
392*860dc821SJeremy L Thompson      write (*, '(A7, I8)', advance='no') ',  nz = ', num_xyz(3)
393*860dc821SJeremy L Thompson    end if
394*860dc821SJeremy L Thompson    write (*, *)
395*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
396*860dc821SJeremy L Thompson    endif
397*860dc821SJeremy L Thompson
398*860dc821SJeremy L Thompson! Build CeedElemRestriction objects describing the mesh and solution discrete representation
399*860dc821SJeremy L Thompson    call buildcartesianrestriction(ceed, fe_dim, num_xyz, mesh_degree, num_comp_x, mesh_size, num_qpts,&
400*860dc821SJeremy L Thompson             &mesh_restriction, ceed_qfunction_none, err)
401*860dc821SJeremy L Thompson    call buildcartesianrestriction(ceed, fe_dim, num_xyz, sol_degree, 1, sol_size, num_qpts, sol_restriction,&
402*860dc821SJeremy L Thompson             &q_data_restriction, err)
403*860dc821SJeremy L Thompson
404*860dc821SJeremy L Thompson    if (test == 0) then
405*860dc821SJeremy L Thompson! LCOV_EXCL_START
406*860dc821SJeremy L Thompson      write (*, *) 'Number of mesh nodes     : ', mesh_size/fe_dim
407*860dc821SJeremy L Thompson      write (*, *) 'Number of solution nodes : ', sol_size
408*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
409*860dc821SJeremy L Thompson    end if
410*860dc821SJeremy L Thompson
411*860dc821SJeremy L Thompson! Create a CeedVector with the mesh coordinates
412*860dc821SJeremy L Thompson! Apply a transformation to the mesh
413*860dc821SJeremy L Thompson    call ceedvectorcreate(ceed, mesh_size, mesh_coords, err)
414*860dc821SJeremy L Thompson    call setcartesianmeshcoords(fe_dim, num_xyz, mesh_degree, mesh_coords, exact_volume, err)
415*860dc821SJeremy L Thompson
416*860dc821SJeremy L Thompson! Context data to be passed to the 'build_mass' QFunction
417*860dc821SJeremy L Thompson    build_ctx_data(1) = fe_dim
418*860dc821SJeremy L Thompson    build_ctx_data(2) = num_comp_x
419*860dc821SJeremy L Thompson    call ceedqfunctioncontextcreate(ceed, build_ctx, err)
420*860dc821SJeremy L Thompson! Note: The context technically only takes arrays of double precision values, but we can pass arrays of ints of the same length
421*860dc821SJeremy L Thompson    offset = 0
422*860dc821SJeremy L Thompson    call ceedqfunctioncontextsetdata(build_ctx, ceed_mem_host, ceed_use_pointer, build_ctx_size, build_ctx_data,&
423*860dc821SJeremy L Thompson             &offset, err)
424*860dc821SJeremy L Thompson
425*860dc821SJeremy L Thompson! Create the QFunction that builds the mass operator (i.e. computes its quadrature data) and set its context data
426*860dc821SJeremy L Thompson    if (gallery == 1) then
427*860dc821SJeremy L Thompson      select case (fe_dim)
428*860dc821SJeremy L Thompson        case (1)
429*860dc821SJeremy L Thompson          call ceedqfunctioncreateinteriorbyname(ceed, 'Mass1DBuild', qf_build, err)
430*860dc821SJeremy L Thompson
431*860dc821SJeremy L Thompson        case (2)
432*860dc821SJeremy L Thompson          call ceedqfunctioncreateinteriorbyname(ceed, 'Mass2DBuild', qf_build, err)
433*860dc821SJeremy L Thompson
434*860dc821SJeremy L Thompson        case (3)
435*860dc821SJeremy L Thompson          call ceedqfunctioncreateinteriorbyname(ceed, 'Mass3DBuild', qf_build, err)
436*860dc821SJeremy L Thompson      end select
437*860dc821SJeremy L Thompson    else
438*860dc821SJeremy L Thompson      call ceedqfunctioncreateinterior(ceed, 1, build_mass,&
439*860dc821SJeremy L Thompson               &SOURCE_DIR&
440*860dc821SJeremy L Thompson               &//'ex1-volume-f-c.h:build_mass'//char(0), qf_build, err)
441*860dc821SJeremy L Thompson      call ceedqfunctionaddinput(qf_build, 'dx', num_comp_x * fe_dim, ceed_eval_grad, err)
442*860dc821SJeremy L Thompson      call ceedqfunctionaddinput(qf_build, 'weights', 1, ceed_eval_weight, err)
443*860dc821SJeremy L Thompson      call ceedqfunctionaddoutput(qf_build, 'qdata', 1, ceed_eval_none, err)
444*860dc821SJeremy L Thompson      call ceedqfunctionsetcontext(qf_build, build_ctx, err)
445*860dc821SJeremy L Thompson    end if
446*860dc821SJeremy L Thompson
447*860dc821SJeremy L Thompson! Create the operator that builds the quadrature data for the mass operator
448*860dc821SJeremy L Thompson    call ceedoperatorcreate(ceed, qf_build, ceed_qfunction_none, ceed_qfunction_none, op_build, err)
449*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_build, 'dx', mesh_restriction, mesh_basis, ceed_vector_active, err)
450*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_build, 'weights', ceed_elemrestriction_none, mesh_basis, ceed_vector_none, err)
451*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_build, 'qdata', q_data_restriction, ceed_basis_none, ceed_vector_active, err)
452*860dc821SJeremy L Thompson
453*860dc821SJeremy L Thompson! Compute the quadrature data for the mass operator
454*860dc821SJeremy L Thompson    num_elem  = 1
455*860dc821SJeremy L Thompson    elem_qpts = num_qpts**fe_dim
456*860dc821SJeremy L Thompson    do i = 1, fe_dim
457*860dc821SJeremy L Thompson      num_elem = num_elem * num_xyz(i)
458*860dc821SJeremy L Thompson    end do
459*860dc821SJeremy L Thompson    call ceedvectorcreate(ceed, num_elem * elem_qpts, q_data, err)
460*860dc821SJeremy L Thompson    call ceedoperatorapply(op_build, mesh_coords, q_data, ceed_request_immediate, err)
461*860dc821SJeremy L Thompson
462*860dc821SJeremy L Thompson! Create the QFunction that defines the action of the mass operator
463*860dc821SJeremy L Thompson    if (gallery == 1) then
464*860dc821SJeremy L Thompson      call ceedqfunctioncreateinteriorbyname(ceed, 'MassApply', qf_apply, err)
465*860dc821SJeremy L Thompson    else
466*860dc821SJeremy L Thompson      call ceedqfunctioncreateinterior(ceed, 1, apply_mass,&
467*860dc821SJeremy L Thompson               &SOURCE_DIR&
468*860dc821SJeremy L Thompson               &//'ex1-volume-f-c.h:apply_mass'//char(0), qf_apply, err)
469*860dc821SJeremy L Thompson      call ceedqfunctionaddinput(qf_apply, 'u', 1, ceed_eval_interp, err)
470*860dc821SJeremy L Thompson      call ceedqfunctionaddinput(qf_apply, 'qdata', 1, ceed_eval_none, err)
471*860dc821SJeremy L Thompson      call ceedqfunctionaddoutput(qf_apply, 'v', 1, ceed_eval_interp, err)
472*860dc821SJeremy L Thompson    end if
473*860dc821SJeremy L Thompson
474*860dc821SJeremy L Thompson! Create the mass operator
475*860dc821SJeremy L Thompson    call ceedoperatorcreate(ceed, qf_apply, ceed_qfunction_none, ceed_qfunction_none, op_apply, err)
476*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_apply, 'u', sol_restriction, sol_basis, ceed_vector_active, err)
477*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_apply, 'qdata', q_data_restriction, ceed_basis_none, q_data, err)
478*860dc821SJeremy L Thompson    call ceedoperatorsetfield(op_apply, 'v', sol_restriction, sol_basis, ceed_vector_active, err)
479*860dc821SJeremy L Thompson
480*860dc821SJeremy L Thompson! Create auxiliary solution-size vectors
481*860dc821SJeremy L Thompson    allocate (u_array(sol_size))
482*860dc821SJeremy L Thompson    allocate (v_array(sol_size))
483*860dc821SJeremy L Thompson
484*860dc821SJeremy L Thompson    call ceedvectorcreate(ceed, sol_size, u, err)
485*860dc821SJeremy L Thompson    offset = 0
486*860dc821SJeremy L Thompson    call ceedvectorsetarray(u, ceed_mem_host, ceed_use_pointer, u_array, offset, err)
487*860dc821SJeremy L Thompson    call ceedvectorcreate(ceed, sol_size, v, err)
488*860dc821SJeremy L Thompson    offset = 0
489*860dc821SJeremy L Thompson    call ceedvectorsetarray(v, ceed_mem_host, ceed_use_pointer, v_array, offset, err)
490*860dc821SJeremy L Thompson
491*860dc821SJeremy L Thompson! Initialize 'u' with ones
492*860dc821SJeremy L Thompson    call ceedvectorsetvalue(u, 1.d0, err)
493*860dc821SJeremy L Thompson
494*860dc821SJeremy L Thompson! Compute the mesh volume using the mass operator: volume = 1^T \cdot M \cdot 1
495*860dc821SJeremy L Thompson    call ceedoperatorapply(op_apply, u, v, ceed_request_immediate, err)
496*860dc821SJeremy L Thompson
497*860dc821SJeremy L Thompson! Benchmark runs
498*860dc821SJeremy L Thompson    if (test /= 1 .AND. benchmark /= 0) then
499*860dc821SJeremy L Thompson! LCOV_EXCL_START
500*860dc821SJeremy L Thompson      write (*, *) ' Executing ', benchmark, ' benchmarking runs...'
501*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
502*860dc821SJeremy L Thompson    end if
503*860dc821SJeremy L Thompson    do i = 1, benchmark
504*860dc821SJeremy L Thompson! LCOV_EXCL_START
505*860dc821SJeremy L Thompson      call ceedoperatorapply(op_apply, u, v, ceed_request_immediate, err)
506*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
507*860dc821SJeremy L Thompson    end do
508*860dc821SJeremy L Thompson
509*860dc821SJeremy L Thompson! Compute and print the sum of the entries of 'v' giving the mesh volume
510*860dc821SJeremy L Thompson    computed_volume = 0.d0
511*860dc821SJeremy L Thompson
512*860dc821SJeremy L Thompson    call ceedvectorgetarrayread(v, ceed_mem_host, v_array, offset, err)
513*860dc821SJeremy L Thompson    do i = 1, sol_size
514*860dc821SJeremy L Thompson      computed_volume = computed_volume + v_array(offset + i)
515*860dc821SJeremy L Thompson    end do
516*860dc821SJeremy L Thompson    call ceedvectorrestorearrayread(v, v_array, offset, err)
517*860dc821SJeremy L Thompson
518*860dc821SJeremy L Thompson    if (test /= 1) then
519*860dc821SJeremy L Thompson! LCOV_EXCL_START
520*860dc821SJeremy L Thompson      write (*, *) ' done.'
521*860dc821SJeremy L Thompson      write (*, *) 'Exact mesh volume    :', exact_volume
522*860dc821SJeremy L Thompson      write (*, *) 'Computed mesh volume :', computed_volume
523*860dc821SJeremy L Thompson      write (*, *) 'Volume error         :', (exact_volume - computed_volume)
524*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
525*860dc821SJeremy L Thompson    else
526*860dc821SJeremy L Thompson      if (fe_dim == 1) then
527*860dc821SJeremy L Thompson        if (abs(exact_volume - computed_volume) > 200.d0 * 1e-15) then
528*860dc821SJeremy L Thompson! LCOV_EXCL_START
529*860dc821SJeremy L Thompson          write (*, *) 'Volume error : ', (exact_volume - computed_volume)
530*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
531*860dc821SJeremy L Thompson        end if
532*860dc821SJeremy L Thompson      else
533*860dc821SJeremy L Thompson        if (abs(exact_volume - computed_volume) > 1e-5) then
534*860dc821SJeremy L Thompson! LCOV_EXCL_START
535*860dc821SJeremy L Thompson          write (*, *) 'Volume error : ', (exact_volume - computed_volume)
536*860dc821SJeremy L Thompson! LCOV_EXCL_STOP
537*860dc821SJeremy L Thompson        end if
538*860dc821SJeremy L Thompson      end if
539*860dc821SJeremy L Thompson    end if
540*860dc821SJeremy L Thompson
541*860dc821SJeremy L Thompson! Free dynamically allocated memory
542*860dc821SJeremy L Thompson    call ceedvectordestroy(mesh_coords, err)
543*860dc821SJeremy L Thompson    call ceedvectordestroy(q_data, err)
544*860dc821SJeremy L Thompson    call ceedvectordestroy(u, err)
545*860dc821SJeremy L Thompson    call ceedvectordestroy(v, err)
546*860dc821SJeremy L Thompson    deallocate (u_array)
547*860dc821SJeremy L Thompson    deallocate (v_array)
548*860dc821SJeremy L Thompson    call ceedbasisdestroy(sol_basis, err)
549*860dc821SJeremy L Thompson    call ceedbasisdestroy(mesh_basis, err)
550*860dc821SJeremy L Thompson    call ceedqfunctioncontextdestroy(build_ctx, err)
551*860dc821SJeremy L Thompson    call ceedqfunctiondestroy(qf_build, err)
552*860dc821SJeremy L Thompson    call ceedqfunctiondestroy(qf_apply, err)
553*860dc821SJeremy L Thompson    call ceedoperatordestroy(op_build, err)
554*860dc821SJeremy L Thompson    call ceedoperatordestroy(op_apply, err)
555*860dc821SJeremy L Thompson    call ceeddestroy(ceed, err)
556*860dc821SJeremy L Thompsonend
557*860dc821SJeremy L Thompson!-----------------------------------------------------------------------
558