xref: /libCEED/examples/fluids/navierstokes.c (revision 1e150236cb75647bab4722af963373560dd08cc4)
1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4ccaff030SJeremy L Thompson //
5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Navier-Stokes
18ccaff030SJeremy L Thompson //
19ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve a
20ccaff030SJeremy L Thompson // Navier-Stokes problem.
21ccaff030SJeremy L Thompson //
22ccaff030SJeremy L Thompson // The code is intentionally "raw", using only low-level communication
23ccaff030SJeremy L Thompson // primitives.
24ccaff030SJeremy L Thompson //
25ccaff030SJeremy L Thompson // Build with:
26ccaff030SJeremy L Thompson //
27ccaff030SJeremy L Thompson //     make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes
28ccaff030SJeremy L Thompson //
29ccaff030SJeremy L Thompson // Sample runs:
30ccaff030SJeremy L Thompson //
3184d34d69SLeila Ghaffari //     ./navierstokes -ceed /cpu/self -problem density_current -degree 1
3284d34d69SLeila Ghaffari //     ./navierstokes -ceed /gpu/occa -problem advection -degree 1
33ccaff030SJeremy L Thompson //
3484d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test explicit -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3
3584d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_stab_none -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha
3684d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_stab_supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -stab supg
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson /// @file
39ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson #include <petscts.h>
44ccaff030SJeremy L Thompson #include <petscdmplex.h>
45ccaff030SJeremy L Thompson #include <ceed.h>
46ccaff030SJeremy L Thompson #include <stdbool.h>
47ccaff030SJeremy L Thompson #include <petscsys.h>
48ccaff030SJeremy L Thompson #include "common.h"
49b0137797SLeila Ghaffari #include "setup-boundary.h"
50ccaff030SJeremy L Thompson #include "advection.h"
51ccaff030SJeremy L Thompson #include "advection2d.h"
52ccaff030SJeremy L Thompson #include "densitycurrent.h"
53ccaff030SJeremy L Thompson 
5484d34d69SLeila Ghaffari #if PETSC_VERSION_LT(3,14,0)
5584d34d69SLeila Ghaffari #  define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i)
5684d34d69SLeila Ghaffari #  define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i)
5784d34d69SLeila Ghaffari #endif
5884d34d69SLeila Ghaffari 
5984d34d69SLeila Ghaffari // MemType Options
6084d34d69SLeila Ghaffari static const char *const memTypes[] = {
6184d34d69SLeila Ghaffari   "host",
6284d34d69SLeila Ghaffari   "device",
6384d34d69SLeila Ghaffari   "memType", "CEED_MEM_", NULL
6484d34d69SLeila Ghaffari };
6584d34d69SLeila Ghaffari 
66ccaff030SJeremy L Thompson // Problem Options
67ccaff030SJeremy L Thompson typedef enum {
68ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
69ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
70ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
71ccaff030SJeremy L Thompson } problemType;
72ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
73ccaff030SJeremy L Thompson   "density_current",
74ccaff030SJeremy L Thompson   "advection",
75ccaff030SJeremy L Thompson   "advection2d",
7684d34d69SLeila Ghaffari   "problemType", "NS_", NULL
77ccaff030SJeremy L Thompson };
78ccaff030SJeremy L Thompson 
79ccaff030SJeremy L Thompson typedef enum {
80ccaff030SJeremy L Thompson   STAB_NONE = 0,
81ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
82ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
83ccaff030SJeremy L Thompson } StabilizationType;
84ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
8584d34d69SLeila Ghaffari   "none",
86ccaff030SJeremy L Thompson   "SU",
87ccaff030SJeremy L Thompson   "SUPG",
88ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
89ccaff030SJeremy L Thompson };
90ccaff030SJeremy L Thompson 
9184d34d69SLeila Ghaffari // Test Options
9284d34d69SLeila Ghaffari typedef enum {
9384d34d69SLeila Ghaffari   TEST_NONE = 0,               // Non test mode
9484d34d69SLeila Ghaffari   TEST_EXPLICIT = 1,           // Explicit test
9584d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab
9684d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab
9784d34d69SLeila Ghaffari } testType;
9884d34d69SLeila Ghaffari static const char *const testTypes[] = {
9984d34d69SLeila Ghaffari   "none",
10084d34d69SLeila Ghaffari   "explicit",
10184d34d69SLeila Ghaffari   "implicit_stab_none",
10284d34d69SLeila Ghaffari   "implicit_stab_supg",
10384d34d69SLeila Ghaffari   "testType", "TEST_", NULL
10484d34d69SLeila Ghaffari };
10584d34d69SLeila Ghaffari 
10684d34d69SLeila Ghaffari // Tests specific data
10784d34d69SLeila Ghaffari typedef struct {
10884d34d69SLeila Ghaffari   PetscScalar testtol;
10984d34d69SLeila Ghaffari   const char *filepath;
11084d34d69SLeila Ghaffari } testData;
11184d34d69SLeila Ghaffari 
11284d34d69SLeila Ghaffari testData testOptions[] = {
11384d34d69SLeila Ghaffari   [TEST_NONE] = {
11484d34d69SLeila Ghaffari     .testtol = 0.,
11584d34d69SLeila Ghaffari     .filepath = NULL
11684d34d69SLeila Ghaffari   },
11784d34d69SLeila Ghaffari   [TEST_EXPLICIT] = {
11884d34d69SLeila Ghaffari     .testtol = 1E-5,
11984d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin"
12084d34d69SLeila Ghaffari   },
12184d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_NONE] = {
12284d34d69SLeila Ghaffari     .testtol = 5E-4,
12384d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin"
12484d34d69SLeila Ghaffari   },
12584d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_SUPG] = {
12684d34d69SLeila Ghaffari     .testtol = 5E-4,
12784d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin"
12884d34d69SLeila Ghaffari   }
12984d34d69SLeila Ghaffari };
13084d34d69SLeila Ghaffari 
131ccaff030SJeremy L Thompson // Problem specific data
132ccaff030SJeremy L Thompson typedef struct {
1338b982baeSLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
134ea6e0f84SLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applySur_rhs,
135ea6e0f84SLeila Ghaffari                     applyVol_ifunction, applySur_ifunction;
136ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
137ccaff030SJeremy L Thompson                        PetscScalar[], void *);
138ea6e0f84SLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
139ea6e0f84SLeila Ghaffari              *applySur_rhs_loc, *applyVol_ifunction_loc, *applySur_ifunction_loc;
140ccaff030SJeremy L Thompson   const bool non_zero_time;
141ccaff030SJeremy L Thompson } problemData;
142ccaff030SJeremy L Thompson 
143ccaff030SJeremy L Thompson problemData problemOptions[] = {
144ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
145ccaff030SJeremy L Thompson     .dim                       = 3,
146ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1478b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
148b0137797SLeila Ghaffari     .setupVol                  = Setup,
149b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
150356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
151356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
152ccaff030SJeremy L Thompson     .ics                       = ICsDC,
153ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
154c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
155c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
1565e197d05SLeila Ghaffari   //.applySur_rhs              = DC_Sur,
1575e197d05SLeila Ghaffari   //.applySur_rhs_loc          = DC_Sur_loc,
158c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
159c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
160ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_DC_Sur,
161ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_DC_Sur_loc,
162ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
16384d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
164ccaff030SJeremy L Thompson   },
165ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
166ccaff030SJeremy L Thompson     .dim                       = 3,
167ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1688b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
169b0137797SLeila Ghaffari     .setupVol                  = Setup,
170b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
171356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
172356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
173ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
174ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
175c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
176c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
17729448992SLeila Ghaffari     .applySur_rhs              = Advection_Sur,
17829448992SLeila Ghaffari     .applySur_rhs_loc          = Advection_Sur_loc,
179c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
180c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
181ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection_Sur,
182ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection_Sur_loc,
183ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
18484d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
185ccaff030SJeremy L Thompson   },
186ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
187ccaff030SJeremy L Thompson     .dim                       = 2,
188ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
1898b982baeSLeila Ghaffari     .qdatasizeSur              = 3,
190c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
191c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
192b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
193b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
194ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
195ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
196c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
197c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
198b0137797SLeila Ghaffari     .applySur_rhs              = Advection2d_Sur,
199b0137797SLeila Ghaffari     .applySur_rhs_loc          = Advection2d_Sur_loc,
200c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
201c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
202ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection2d_Sur,
203ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection2d_Sur_loc,
204ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
20584d34d69SLeila Ghaffari     .non_zero_time             = PETSC_TRUE,
206ccaff030SJeremy L Thompson   },
207ccaff030SJeremy L Thompson };
208ccaff030SJeremy L Thompson 
209ccaff030SJeremy L Thompson // PETSc user data
210ccaff030SJeremy L Thompson typedef struct User_ *User;
211ccaff030SJeremy L Thompson typedef struct Units_ *Units;
212ccaff030SJeremy L Thompson 
213ccaff030SJeremy L Thompson struct User_ {
214ccaff030SJeremy L Thompson   MPI_Comm comm;
215ccaff030SJeremy L Thompson   PetscInt outputfreq;
216ccaff030SJeremy L Thompson   DM dm;
217ccaff030SJeremy L Thompson   DM dmviz;
218ccaff030SJeremy L Thompson   Mat interpviz;
219ccaff030SJeremy L Thompson   Ceed ceed;
220ccaff030SJeremy L Thompson   Units units;
221ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
222*1e150236SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction;
223ccaff030SJeremy L Thompson   Vec M;
224ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
225ccaff030SJeremy L Thompson   PetscInt contsteps;
226ccaff030SJeremy L Thompson };
227ccaff030SJeremy L Thompson 
228ccaff030SJeremy L Thompson struct Units_ {
229ccaff030SJeremy L Thompson   // fundamental units
230ccaff030SJeremy L Thompson   PetscScalar meter;
231ccaff030SJeremy L Thompson   PetscScalar kilogram;
232ccaff030SJeremy L Thompson   PetscScalar second;
233ccaff030SJeremy L Thompson   PetscScalar Kelvin;
234ccaff030SJeremy L Thompson   // derived units
235ccaff030SJeremy L Thompson   PetscScalar Pascal;
236ccaff030SJeremy L Thompson   PetscScalar JperkgK;
237ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
238ccaff030SJeremy L Thompson   PetscScalar WpermK;
239ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
240ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
241ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
242ccaff030SJeremy L Thompson };
243ccaff030SJeremy L Thompson 
244ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
245ccaff030SJeremy L Thompson struct SimpleBC_ {
246cfa64770SLeila Ghaffari   PetscInt nwall, nslip[3], noutflow;
24784d34d69SLeila Ghaffari   PetscInt walls[6], slips[3][6], outflow[6];
24884d34d69SLeila Ghaffari   PetscBool userbc;
249ccaff030SJeremy L Thompson };
250ccaff030SJeremy L Thompson 
251ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
252ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
253ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
254ccaff030SJeremy L Thompson }
255ccaff030SJeremy L Thompson 
256ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
257ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
25884d34d69SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
259ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
260ccaff030SJeremy L Thompson 
261ccaff030SJeremy L Thompson   PetscSection   section;
2620c6c0b13SLeila Ghaffari   PetscInt       p, Nelem, Ndof, *erestrict, eoffset, nfields, dim,
2630c6c0b13SLeila Ghaffari                  depth;
2640c6c0b13SLeila Ghaffari   DMLabel        depthLabel;
2650c6c0b13SLeila Ghaffari   IS             depthIS, iterIS;
26684d34d69SLeila Ghaffari   Vec            Uloc;
2670c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
268ccaff030SJeremy L Thompson   PetscErrorCode ierr;
269ccaff030SJeremy L Thompson 
270ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
271ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
272da51bdd9SLeila Ghaffari   dim -= height;
273ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
274ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
275ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
276ccaff030SJeremy L Thompson   fieldoff[0] = 0;
277ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
278ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
279ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
280ccaff030SJeremy L Thompson   }
281ccaff030SJeremy L Thompson 
2820c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2830c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2840c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2850c6c0b13SLeila Ghaffari   if (domainLabel) {
2860c6c0b13SLeila Ghaffari     IS domainIS;
2870c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2880c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2890c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2900c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2910c6c0b13SLeila Ghaffari   } else {
2920c6c0b13SLeila Ghaffari     iterIS = depthIS;
2930c6c0b13SLeila Ghaffari   }
2940c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2950c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
296ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2970c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2980c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
299ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
30084d34d69SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
30184d34d69SLeila Ghaffari                                    &numindices, &indices, NULL, NULL);
30284d34d69SLeila Ghaffari     CHKERRQ(ierr);
30384d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
30484d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
30584d34d69SLeila Ghaffari           c);
306ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
307ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
308ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
309ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
310ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
311ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
312ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
313ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
314ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
315ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
316ccaff030SJeremy L Thompson                      c, i, f, j);
317ccaff030SJeremy L Thompson         }
318ccaff030SJeremy L Thompson       }
319ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
320ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
3216f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
322ccaff030SJeremy L Thompson     }
32384d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
32484d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
32584d34d69SLeila Ghaffari     CHKERRQ(ierr);
326ccaff030SJeremy L Thompson   }
3270c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3280c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3290c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
330ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3310c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3320c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3330c6c0b13SLeila Ghaffari 
334ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
335ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
336ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3376f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3386f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3396f55dfd5Svaleriabarra                             Erestrict);
340ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
341ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
342ccaff030SJeremy L Thompson }
343ccaff030SJeremy L Thompson 
344c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
345*1e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
346*1e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
347*1e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
348*1e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
349c96c872fSLeila Ghaffari 
350c96c872fSLeila Ghaffari   DM dmcoord;
351*1e150236SLeila Ghaffari   CeedInt dim, localNelem;
352*1e150236SLeila Ghaffari   CeedInt Qdim;
353c96c872fSLeila Ghaffari   PetscErrorCode ierr;
354c96c872fSLeila Ghaffari 
355c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
356*1e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
357*1e150236SLeila Ghaffari   dim -= height;
358*1e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
359c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
360c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
361c96c872fSLeila Ghaffari   CHKERRQ(ierr);
362c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
363c96c872fSLeila Ghaffari   CHKERRQ(ierr);
364c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
365c96c872fSLeila Ghaffari   CHKERRQ(ierr);
366c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
367c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
368c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
369c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
370c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
371c96c872fSLeila Ghaffari }
372c96c872fSLeila Ghaffari 
373*1e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
374*1e150236SLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, CeedOperator op_applyVol,
375ca3ac6ddSLeila Ghaffari     CeedQFunction qf_applySur, CeedQFunction qf_setupSur, CeedInt height, PetscInt nFace,
376ca3ac6ddSLeila Ghaffari     PetscInt value[6], CeedInt numP_Sur, CeedInt numQ_Sur, CeedInt qdatasizeSur, CeedInt NqptsSur,
377*1e150236SLeila Ghaffari     CeedBasis basisxSur, CeedBasis basisqSur, CeedOperator *op_apply) {
378ca3ac6ddSLeila Ghaffari 
379ca3ac6ddSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
380ca3ac6ddSLeila Ghaffari   PetscInt dim, localNelemSur[6];
381*1e150236SLeila Ghaffari   Vec Xloc;
382*1e150236SLeila Ghaffari   CeedVector xcorners, qdataSur[6];
383ca3ac6ddSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
384ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
385*1e150236SLeila Ghaffari   PetscScalar *x;
386*1e150236SLeila Ghaffari   PetscInt lsize;
387ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
388ca3ac6ddSLeila Ghaffari 
389ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
390ca3ac6ddSLeila Ghaffari   // Composite Operaters
391ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
392*1e150236SLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_applyVol); // Apply a Sub-Operator for the volume
393ca3ac6ddSLeila Ghaffari 
394ca3ac6ddSLeila Ghaffari   if (nFace) {
395*1e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
396*1e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
397*1e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
398*1e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
399*1e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
400*1e150236SLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
401ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
402*1e150236SLeila Ghaffari     // Create CEED Operator for each In/OutFlow faces
403ca3ac6ddSLeila Ghaffari     for(CeedInt i=0; i<nFace; i++) {
404*1e150236SLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, value[i], numP_Sur,
405ca3ac6ddSLeila Ghaffari                                      numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
406ca3ac6ddSLeila Ghaffari                                      &restrictqdiSur[i]); CHKERRQ(ierr);
407*1e150236SLeila Ghaffari 
408*1e150236SLeila Ghaffari       // Create the CEED vectors that will be needed in boundary setup
409ca3ac6ddSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
410ca3ac6ddSLeila Ghaffari       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
411ca3ac6ddSLeila Ghaffari 
412*1e150236SLeila Ghaffari       // Create the operator that builds the quadrature data for the In/OutFlow operator
413ca3ac6ddSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
414ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, CEED_VECTOR_ACTIVE);
415ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
416ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
417ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
418ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
419*1e150236SLeila Ghaffari 
420*1e150236SLeila Ghaffari       // Create In/OutFlow operator
421ca3ac6ddSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
422ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
423ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
424ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
425ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, xcorners);
426ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
427ca3ac6ddSLeila Ghaffari 
428*1e150236SLeila Ghaffari       // Apply CEED operator for boundary setup
429*1e150236SLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i], CEED_REQUEST_IMMEDIATE);
430*1e150236SLeila Ghaffari 
431*1e150236SLeila Ghaffari       // Apply Sub-Operator for In/OutFlow BCs
432ca3ac6ddSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
433ca3ac6ddSLeila Ghaffari     }
434*1e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
435ca3ac6ddSLeila Ghaffari   }
436ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
437ca3ac6ddSLeila Ghaffari }
438ca3ac6ddSLeila Ghaffari 
439ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
440ccaff030SJeremy L Thompson   PetscErrorCode ierr;
441ccaff030SJeremy L Thompson   PetscInt m;
442ccaff030SJeremy L Thompson 
443ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
444ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
445ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
446ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
447ccaff030SJeremy L Thompson }
448ccaff030SJeremy L Thompson 
449ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
450ccaff030SJeremy L Thompson   PetscErrorCode ierr;
451ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
452ccaff030SJeremy L Thompson   PetscScalar *a;
453ccaff030SJeremy L Thompson 
454ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
455ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
456ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
457ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
45884d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
45984d34d69SLeila Ghaffari                                   mpetsc, mceed);
460ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
461ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
462ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
463ccaff030SJeremy L Thompson }
464ccaff030SJeremy L Thompson 
465ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
466ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
467ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
468ccaff030SJeremy L Thompson   PetscErrorCode ierr;
469ccaff030SJeremy L Thompson   Vec Qbc;
470ccaff030SJeremy L Thompson 
471ccaff030SJeremy L Thompson   PetscFunctionBegin;
472ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
473ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
474ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
475ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
476ccaff030SJeremy L Thompson }
477ccaff030SJeremy L Thompson 
478ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
479ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
480ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
481ccaff030SJeremy L Thompson   PetscErrorCode ierr;
482ccaff030SJeremy L Thompson   User user = *(User *)userData;
483ccaff030SJeremy L Thompson   PetscScalar *q, *g;
484ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
485ccaff030SJeremy L Thompson 
486ccaff030SJeremy L Thompson   // Global-to-local
487ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
488ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
489ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
490ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
491ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
493ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
494ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
495ccaff030SJeremy L Thompson 
496ccaff030SJeremy L Thompson   // Ceed Vectors
497ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
499ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
500ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
501ccaff030SJeremy L Thompson 
502ccaff030SJeremy L Thompson   // Apply CEED operator
503ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
504ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
505ccaff030SJeremy L Thompson 
506ccaff030SJeremy L Thompson   // Restore vectors
507ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
508ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
509ccaff030SJeremy L Thompson 
510ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
511ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
512ccaff030SJeremy L Thompson 
513ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
514ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
515ccaff030SJeremy L Thompson   CHKERRQ(ierr);
516ccaff030SJeremy L Thompson 
517ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
518ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
520ccaff030SJeremy L Thompson }
521ccaff030SJeremy L Thompson 
522ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
523ccaff030SJeremy L Thompson                                    void *userData) {
524ccaff030SJeremy L Thompson   PetscErrorCode ierr;
525ccaff030SJeremy L Thompson   User user = *(User *)userData;
526ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
527ccaff030SJeremy L Thompson   PetscScalar *g;
528ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
529ccaff030SJeremy L Thompson 
530ccaff030SJeremy L Thompson   // Global-to-local
531ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
532ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
536ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
538ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
539ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson 
543ccaff030SJeremy L Thompson   // Ceed Vectors
544ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
545ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
546ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
548ccaff030SJeremy L Thompson                      (PetscScalar *)q);
549ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
550ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
551ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson   // Apply CEED operator
554ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
555ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
556ccaff030SJeremy L Thompson 
557ccaff030SJeremy L Thompson   // Restore vectors
558ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
560ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson 
562ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
563ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson 
565ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
567ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
568ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
569ccaff030SJeremy L Thompson }
570ccaff030SJeremy L Thompson 
571ccaff030SJeremy L Thompson // User provided TS Monitor
572ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
573ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
574ccaff030SJeremy L Thompson   User user = ctx;
575ccaff030SJeremy L Thompson   Vec Qloc;
576ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
577ccaff030SJeremy L Thompson   PetscViewer viewer;
578ccaff030SJeremy L Thompson   PetscErrorCode ierr;
579ccaff030SJeremy L Thompson 
580ccaff030SJeremy L Thompson   // Set up output
581ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
582ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
583ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
584ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
585ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
586ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
589ccaff030SJeremy L Thompson 
590ccaff030SJeremy L Thompson   // Output
591ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
592ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
593ccaff030SJeremy L Thompson   CHKERRQ(ierr);
594ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
595ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
596ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
5979d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
598ccaff030SJeremy L Thompson   if (user->dmviz) {
599ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
600ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
601ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
602ccaff030SJeremy L Thompson 
603ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
604ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
605ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
606ccaff030SJeremy L Thompson     CHKERRQ(ierr);
607ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
609ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
610ccaff030SJeremy L Thompson     CHKERRQ(ierr);
611ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
612ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
613ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
614ccaff030SJeremy L Thompson     CHKERRQ(ierr);
615ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
616ccaff030SJeremy L Thompson                               filepath_refined,
617ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
618ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
619ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
620ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
621ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
622ccaff030SJeremy L Thompson   }
623ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
624ccaff030SJeremy L Thompson 
625ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
626ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
627ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
628ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
629ccaff030SJeremy L Thompson   CHKERRQ(ierr);
630ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
631ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
632ccaff030SJeremy L Thompson 
633ccaff030SJeremy L Thompson   // Save time stamp
634ccaff030SJeremy L Thompson   // Dimensionalize time back
635ccaff030SJeremy L Thompson   time /= user->units->second;
636ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
637ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
638ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
639ccaff030SJeremy L Thompson   CHKERRQ(ierr);
640ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
641ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
642ccaff030SJeremy L Thompson   #else
643ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
644ccaff030SJeremy L Thompson   #endif
645ccaff030SJeremy L Thompson   CHKERRQ(ierr);
646ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
647ccaff030SJeremy L Thompson 
648ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
649ccaff030SJeremy L Thompson }
650ccaff030SJeremy L Thompson 
65184d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
652ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
653ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
654ccaff030SJeremy L Thompson   PetscErrorCode ierr;
655ccaff030SJeremy L Thompson   CeedVector multlvec;
656ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
657ccaff030SJeremy L Thompson 
658ccaff030SJeremy L Thompson   ctxSetup->time = time;
659ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
660ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
661ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
662ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
663ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
664ccaff030SJeremy L Thompson 
665ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
666ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
667ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
668ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
669ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
670ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
671ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
672ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
673ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
674ccaff030SJeremy L Thompson   CHKERRQ(ierr);
675ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
676ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
677ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
678ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
679ccaff030SJeremy L Thompson 
680ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
681ccaff030SJeremy L Thompson }
682ccaff030SJeremy L Thompson 
683ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
684ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
685ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
686ccaff030SJeremy L Thompson   PetscErrorCode ierr;
687ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
688ccaff030SJeremy L Thompson   CeedOperator op_mass;
689ccaff030SJeremy L Thompson   CeedVector mceed;
690ccaff030SJeremy L Thompson   Vec Mloc;
691ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
692ccaff030SJeremy L Thompson 
693ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
694ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
695ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
696ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
697ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
698ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
699ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
700ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
701ccaff030SJeremy L Thompson 
702ccaff030SJeremy L Thompson   // Create the mass operator
703ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
704ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
705ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
706ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
707ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
708ccaff030SJeremy L Thompson 
709ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
710ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
711ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
712ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
713ccaff030SJeremy L Thompson 
714ccaff030SJeremy L Thompson   {
715ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
716ccaff030SJeremy L Thompson     CeedVector onesvec;
717ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
718ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
719ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
720ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
721ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
722ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
723ccaff030SJeremy L Thompson   }
724ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
725ccaff030SJeremy L Thompson 
726ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
727ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
728ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
729ccaff030SJeremy L Thompson 
730ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
731ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
732ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
733ccaff030SJeremy L Thompson }
734ccaff030SJeremy L Thompson 
73584d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
736ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
737ccaff030SJeremy L Thompson   PetscErrorCode ierr;
738ccaff030SJeremy L Thompson 
739ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
740ccaff030SJeremy L Thompson   {
741ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
742ccaff030SJeremy L Thompson     PetscFE fe;
743ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
744ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
745ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
74632ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
747ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
748ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
749ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
75007af6069Svaleriabarra     {
75107af6069Svaleriabarra       PetscInt comps[1] = {1};
75207af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
75307af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
75407af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
75507af6069Svaleriabarra       comps[0] = 2;
75607af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
75707af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
75807af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
75907af6069Svaleriabarra       comps[0] = 3;
76007af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
76107af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
76207af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
76307af6069Svaleriabarra     }
76484d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
76584d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
76684d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
76784d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
76884d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
76984d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
77084d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
77184d34d69SLeila Ghaffari 
77284d34d69SLeila Ghaffari           }
77384d34d69SLeila Ghaffari         }
77484d34d69SLeila Ghaffari       }
77584d34d69SLeila Ghaffari     }
77684d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
77784d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
77884d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
77984d34d69SLeila Ghaffari     {
78084d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
78184d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
78284d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
78384d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
78484d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
78584d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
78684d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
78784d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
78884d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
78984d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
79084d34d69SLeila Ghaffari       } else
79184d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
79284d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
79384d34d69SLeila Ghaffari     }
794ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
795ccaff030SJeremy L Thompson     CHKERRQ(ierr);
796ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
797ccaff030SJeremy L Thompson   }
798ccaff030SJeremy L Thompson   {
799ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
800ccaff030SJeremy L Thompson     PetscSection section;
801ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
802ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
803ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
804ccaff030SJeremy L Thompson     CHKERRQ(ierr);
805ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
806ccaff030SJeremy L Thompson     CHKERRQ(ierr);
807ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
808ccaff030SJeremy L Thompson     CHKERRQ(ierr);
809ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
810ccaff030SJeremy L Thompson     CHKERRQ(ierr);
811ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
812ccaff030SJeremy L Thompson     CHKERRQ(ierr);
813ccaff030SJeremy L Thompson   }
814ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
815ccaff030SJeremy L Thompson }
816ccaff030SJeremy L Thompson 
817ccaff030SJeremy L Thompson int main(int argc, char **argv) {
818ccaff030SJeremy L Thompson   PetscInt ierr;
819ccaff030SJeremy L Thompson   MPI_Comm comm;
82084d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
821ccaff030SJeremy L Thompson   Mat interpviz;
822ccaff030SJeremy L Thompson   TS ts;
823ccaff030SJeremy L Thompson   TSAdapt adapt;
824ccaff030SJeremy L Thompson   User user;
825ccaff030SJeremy L Thompson   Units units;
826ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
82784d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
828ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
829ccaff030SJeremy L Thompson   PetscMPIInt rank;
830ccaff030SJeremy L Thompson   PetscScalar ftime;
831ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
832ccaff030SJeremy L Thompson   Ceed ceed;
83384d34d69SLeila Ghaffari   CeedInt numP, numQ;
834cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
83584d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
83684d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
837cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
838cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
839ccaff030SJeremy L Thompson   CeedScalar Rd;
84084d34d69SLeila Ghaffari   CeedMemType memtyperequested;
841ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
842ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
843ccaff030SJeremy L Thompson   problemType problemChoice;
844ccaff030SJeremy L Thompson   problemData *problem = NULL;
845ccaff030SJeremy L Thompson   StabilizationType stab;
84684d34d69SLeila Ghaffari   testType testChoice;
84784d34d69SLeila Ghaffari   testData *test = NULL;
84884d34d69SLeila Ghaffari   PetscBool implicit;
849cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
850ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
85184d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
85284d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
853ccaff030SJeremy L Thompson   };
854ccaff030SJeremy L Thompson   double start, cpu_time_used;
85584d34d69SLeila Ghaffari   // Check PETSc CUDA support
85684d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
85784d34d69SLeila Ghaffari   // *INDENT-OFF*
85884d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
85984d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
86084d34d69SLeila Ghaffari   #else
86184d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
86284d34d69SLeila Ghaffari   #endif
86384d34d69SLeila Ghaffari   // *INDENT-ON*
864ccaff030SJeremy L Thompson 
865ccaff030SJeremy L Thompson   // Create the libCEED contexts
866ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
867ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
868ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
869ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
870ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
871ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
872ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
873ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
874ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
875ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
876ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
877ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
878ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
879ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
880ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
881ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
882ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
883ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
884ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
885ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
886ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
887ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
888ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
889ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
890ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
891ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
89284d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
89384d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
894ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
895ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
896ccaff030SJeremy L Thompson 
897ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
898ccaff030SJeremy L Thompson   if (ierr) return ierr;
899ccaff030SJeremy L Thompson 
900ccaff030SJeremy L Thompson   // Allocate PETSc context
901ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
902ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
903ccaff030SJeremy L Thompson 
904ccaff030SJeremy L Thompson   // Parse command line options
905ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
906ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
907ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
908ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
909ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
910ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
91184d34d69SLeila Ghaffari   testChoice = TEST_NONE;
91284d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
91384d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
91484d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
91584d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
91684d34d69SLeila Ghaffari   test = &testOptions[testChoice];
917ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
918ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
919ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
920ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
921ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
922ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
923ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
924ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
925ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
926ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
927ccaff030SJeremy L Thompson   CHKERRQ(ierr);
92884d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
92984d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
93084d34d69SLeila Ghaffari     CHKERRQ(ierr);
93184d34d69SLeila Ghaffari   }
932ccaff030SJeremy L Thompson   {
9337573aee6SJed Brown     PetscInt len;
9347573aee6SJed Brown     PetscBool flg;
9357573aee6SJed Brown     ierr = PetscOptionsIntArray("-bc_outflow",
9367573aee6SJed Brown                               "Use outflow boundary conditions on this list of faces",
9377573aee6SJed Brown                               NULL, bc.outflow,
9387573aee6SJed Brown                               (len = sizeof(bc.outflow) / sizeof(bc.outflow[0]),
9397573aee6SJed Brown                               &len), &flg); CHKERRQ(ierr);
9407573aee6SJed Brown     if (flg) {
9417573aee6SJed Brown       bc.noutflow = len;
9427573aee6SJed Brown       // Using outflow boundaries disables automatic wall/slip boundaries (they must be set explicitly)
9437573aee6SJed Brown       bc.nwall = 0;
9447573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
9457573aee6SJed Brown     }
946ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
947ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
948ccaff030SJeremy L Thompson                                 NULL, bc.walls,
949ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
950ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
9517573aee6SJed Brown     if (flg) {
9527573aee6SJed Brown       bc.nwall = len;
9537573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
9547573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
9557573aee6SJed Brown     }
956ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
957ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
958ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
959ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
960ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
961ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
962ccaff030SJeremy L Thompson                                    &len), &flg);
963ccaff030SJeremy L Thompson       CHKERRQ(ierr);
96484d34d69SLeila Ghaffari       if (flg) {
96584d34d69SLeila Ghaffari         bc.nslip[j] = len;
96684d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
96784d34d69SLeila Ghaffari       }
968ccaff030SJeremy L Thompson     }
969ccaff030SJeremy L Thompson   }
970cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
971cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
972cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
973ccaff030SJeremy L Thompson   CHKERRQ(ierr);
974ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
975ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
976ccaff030SJeremy L Thompson   meter = fabs(meter);
977ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
978ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
979ccaff030SJeremy L Thompson   second = fabs(second);
980ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
981ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
982ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
983ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
984ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
985ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
986ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
987ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
988ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
989ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
990ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
991ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
992ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
993ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
994ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
995ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
996ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
997ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
998ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
999ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1000ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1002ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1003ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1004ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1005ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1006ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1007ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1008ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1009ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1010ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
101184d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
101284d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
101384d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
101484d34d69SLeila Ghaffari     CHKERRQ(ierr);
101584d34d69SLeila Ghaffari   }
1016ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1017ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1018ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1019ccaff030SJeremy L Thompson   CHKERRQ(ierr);
102084d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
102184d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
102284d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
102384d34d69SLeila Ghaffari     CHKERRQ(ierr);
102484d34d69SLeila Ghaffari   }
1025ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1026ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1027ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1028ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1029ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1030ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1031ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1032ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1033ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1034ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1035ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1036ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1037ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1038ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
1039ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
1040ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1041ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1042ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1043ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1044ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1045ccaff030SJeremy L Thompson   n = problem->dim;
1046ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1047ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1048ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1049ccaff030SJeremy L Thompson   {
1050ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1051ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1052ccaff030SJeremy L Thompson     if (norm > 0) {
1053ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1054ccaff030SJeremy L Thompson     }
1055ccaff030SJeremy L Thompson   }
1056ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1057ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1058ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1059ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1060ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
106184d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
106284d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
106384d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
106484d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
106584d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1066ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1067ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1068ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
106984d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
107084d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
107184d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
107284d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
107384d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
107484d34d69SLeila Ghaffari   CHKERRQ(ierr);
1075ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1076ccaff030SJeremy L Thompson 
1077ccaff030SJeremy L Thompson   // Define derived units
1078ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1079ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1080ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1081ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1082ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1083ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1084ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
1085ccaff030SJeremy L Thompson 
1086ccaff030SJeremy L Thompson   // Scale variables to desired units
1087ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1088ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1089ccaff030SJeremy L Thompson   P0 *= Pascal;
1090ccaff030SJeremy L Thompson   N *= (1./second);
1091ccaff030SJeremy L Thompson   cv *= JperkgK;
1092ccaff030SJeremy L Thompson   cp *= JperkgK;
1093ccaff030SJeremy L Thompson   Rd = cp - cv;
1094ccaff030SJeremy L Thompson   g *= mpersquareds;
1095ccaff030SJeremy L Thompson   mu *= Pascal * second;
1096ccaff030SJeremy L Thompson   k *= WpermK;
1097ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1098ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1099ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1100ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1101ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1102ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1103ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1104ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1105ccaff030SJeremy L Thompson 
1106ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1107cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1108ccaff030SJeremy L Thompson   // Set up the libCEED context
1109ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1110ccaff030SJeremy L Thompson     .theta0 = theta0,
1111ccaff030SJeremy L Thompson     .thetaC = thetaC,
1112ccaff030SJeremy L Thompson     .P0 = P0,
1113ccaff030SJeremy L Thompson     .N = N,
1114ccaff030SJeremy L Thompson     .cv = cv,
1115ccaff030SJeremy L Thompson     .cp = cp,
1116ccaff030SJeremy L Thompson     .Rd = Rd,
1117ccaff030SJeremy L Thompson     .g = g,
1118ccaff030SJeremy L Thompson     .rc = rc,
1119ccaff030SJeremy L Thompson     .lx = lx,
1120ccaff030SJeremy L Thompson     .ly = ly,
1121ccaff030SJeremy L Thompson     .lz = lz,
1122ccaff030SJeremy L Thompson     .center[0] = center[0],
1123ccaff030SJeremy L Thompson     .center[1] = center[1],
1124ccaff030SJeremy L Thompson     .center[2] = center[2],
1125ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1126ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1127ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
1128ccaff030SJeremy L Thompson     .time = 0,
1129ccaff030SJeremy L Thompson   };
1130ccaff030SJeremy L Thompson 
113184d34d69SLeila Ghaffari   // Create the mesh
1132ccaff030SJeremy L Thompson   {
1133ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1134ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
113584d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1136ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1137ccaff030SJeremy L Thompson   }
113884d34d69SLeila Ghaffari 
113984d34d69SLeila Ghaffari   // Distribute the mesh over processes
114084d34d69SLeila Ghaffari   {
1141ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1142ccaff030SJeremy L Thompson     PetscPartitioner part;
1143ccaff030SJeremy L Thompson 
1144ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1145ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1146ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1147ccaff030SJeremy L Thompson     if (dmDist) {
1148ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1149ccaff030SJeremy L Thompson       dm  = dmDist;
1150ccaff030SJeremy L Thompson     }
1151ccaff030SJeremy L Thompson   }
1152ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1153ccaff030SJeremy L Thompson 
115484d34d69SLeila Ghaffari   // Setup DM
1155ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1156ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
115784d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
115884d34d69SLeila Ghaffari 
115984d34d69SLeila Ghaffari   // Refine DM for high-order viz
1160ccaff030SJeremy L Thompson   dmviz = NULL;
1161ccaff030SJeremy L Thompson   interpviz = NULL;
1162ccaff030SJeremy L Thompson   if (viz_refine) {
1163ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1164ff6701fcSJed Brown 
1165ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1166ff6701fcSJed Brown     dmhierarchy[0] = dm;
116784d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1168ff6701fcSJed Brown       Mat interp_next;
1169ff6701fcSJed Brown 
1170ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1171ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1172ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1173ff6701fcSJed Brown       d = (d + 1) / 2;
1174ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1175ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1176ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1177ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1178ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1179ff6701fcSJed Brown       else {
1180ff6701fcSJed Brown         Mat C;
1181ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1182ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1183ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1184ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1185ff6701fcSJed Brown         interpviz = C;
1186ff6701fcSJed Brown       }
1187ff6701fcSJed Brown     }
1188cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1189ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1190cb3e2689Svaleriabarra     }
1191ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1192ccaff030SJeremy L Thompson   }
1193ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1194ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1195ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1196ccaff030SJeremy L Thompson   lnodes /= ncompq;
1197ccaff030SJeremy L Thompson 
119884d34d69SLeila Ghaffari   // Initialize CEED
119984d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
120084d34d69SLeila Ghaffari   // Set memtype
120184d34d69SLeila Ghaffari   CeedMemType memtypebackend;
120284d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
120384d34d69SLeila Ghaffari   // Check memtype compatibility
120484d34d69SLeila Ghaffari   if (!setmemtyperequest)
120584d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
120684d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
120784d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
120884d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
120984d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
121084d34d69SLeila Ghaffari 
121184d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
121284d34d69SLeila Ghaffari   numP = degree + 1;
121384d34d69SLeila Ghaffari   numQ = numP + qextra;
121484d34d69SLeila Ghaffari 
121584d34d69SLeila Ghaffari     // Print summary
121684d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1217ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1218ccaff030SJeremy L Thompson     int comm_size;
1219ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1220ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1221ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
122284d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1223ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1224ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1225ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
122684d34d69SLeila Ghaffari     const char *usedresource;
122784d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1228ccaff030SJeremy L Thompson 
122984d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
123084d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
123184d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
123284d34d69SLeila Ghaffari                        "  Problem:\n"
123384d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
123484d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
123584d34d69SLeila Ghaffari                        "  PETSc:\n"
123684d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
123784d34d69SLeila Ghaffari                        "  libCEED:\n"
123884d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
123984d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
124084d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
124184d34d69SLeila Ghaffari                        "  Mesh:\n"
124284d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
124384d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
124484d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
124584d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
124684d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
124784d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
124884d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
124984d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
125084d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
125184d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
125284d34d69SLeila Ghaffari                        (setmemtyperequest) ?
125384d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
125484d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
125584d34d69SLeila Ghaffari     CHKERRQ(ierr);
12560c6c0b13SLeila Ghaffari   }
12570c6c0b13SLeila Ghaffari 
1258ccaff030SJeremy L Thompson   // Set up global mass vector
1259ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1260ccaff030SJeremy L Thompson 
126184d34d69SLeila Ghaffari   // Set up libCEED
1262ccaff030SJeremy L Thompson   // CEED Bases
1263ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
126484d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
126584d34d69SLeila Ghaffari                                   &basisq);
126684d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
126784d34d69SLeila Ghaffari                                   &basisx);
126884d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
126984d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1270ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1271ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1272ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1273ccaff030SJeremy L Thompson 
1274ccaff030SJeremy L Thompson   // CEED Restrictions
1275*1e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
127684d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
127784d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1278ccaff030SJeremy L Thompson 
1279ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1280ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1281ccaff030SJeremy L Thompson 
1282ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1283bd910870SLeila Ghaffari   CeedInt NqptsVol;
128484d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
128584d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
12868b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
128784d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1288ccaff030SJeremy L Thompson 
1289ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1290ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1291ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1292ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1293ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
12948b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1295ccaff030SJeremy L Thompson 
1296ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1297ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1298ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1299ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1300ccaff030SJeremy L Thompson 
1301ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1302ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1303ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1304ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1305ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1306ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13078b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1308ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1309ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1310ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1311ccaff030SJeremy L Thompson   }
1312ccaff030SJeremy L Thompson 
1313ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1314ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1315ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1316ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1317ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1318ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1319ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
13208b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1321ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1322ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1323ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1324ccaff030SJeremy L Thompson   }
1325ccaff030SJeremy L Thompson 
1326ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1327ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
132884d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1329ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
133084d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
133184d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1332ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1333ccaff030SJeremy L Thompson 
1334ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1335ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
133684d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
133784d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1338ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1339ccaff030SJeremy L Thompson 
134084d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
134184d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
134284d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1343ccaff030SJeremy L Thompson 
1344ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1345ccaff030SJeremy L Thompson     CeedOperator op;
1346ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
134784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
134884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
134984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
13508b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
135184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
135284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
135384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1354d3630711SJed Brown     user->op_rhs_vol = op;
1355ccaff030SJeremy L Thompson   }
1356ccaff030SJeremy L Thompson 
1357ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1358ccaff030SJeremy L Thompson     CeedOperator op;
1359ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
136084d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
136184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
136284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
136384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
13648b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
136584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
136684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
136784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1368d3630711SJed Brown     user->op_ifunction_vol = op;
1369ccaff030SJeremy L Thompson   }
1370ccaff030SJeremy L Thompson 
1371cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1372cfa64770SLeila Ghaffari   // Outflow Boundary Condition
13736a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
13746a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
13756a0edaf9SLeila Ghaffari   CeedInt height = 1;
13766a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
1377*1e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
1378*1e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1379cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1380cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1381cfa64770SLeila Ghaffari   CeedInt NqptsSur;
1382cfa64770SLeila Ghaffari   CeedQFunction qf_setupSur, qf_rhsSur, qf_ifunctionSur;
1383cfa64770SLeila Ghaffari 
1384cfa64770SLeila Ghaffari   // CEED bases for the boundaries
13856a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
13866a0edaf9SLeila Ghaffari                                   &basisqSur);
13876a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
13886a0edaf9SLeila Ghaffari                                   &basisxSur);
13896a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
13906a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
13916a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
13926a0edaf9SLeila Ghaffari 
1393cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
13946a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
13956a0edaf9SLeila Ghaffari                               &qf_setupSur);
13966a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
13976a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
13986a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
13996a0edaf9SLeila Ghaffari 
14006a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
1401cfa64770SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator on the Surface
14026a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
14036a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
14046a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
14056a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14066a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
14076a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
14086a0edaf9SLeila Ghaffari   }
14096a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
14106a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
14116a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
14126a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
14136a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
14146a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14156a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
14166a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
14176a0edaf9SLeila Ghaffari   }
1418cfa64770SLeila Ghaffari   // Create CEED Operator for each face
1419*1e150236SLeila Ghaffari   if (qf_rhsSur) {
1420*1e150236SLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, user->op_rhs_vol, qf_rhsSur, qf_setupSur, height,
1421ca3ac6ddSLeila Ghaffari                                   bc.noutflow, bc.outflow, numP_Sur, numQ_Sur, qdatasizeSur,
1422*1e150236SLeila Ghaffari                                   NqptsSur, basisxSur, basisqSur, &user->op_rhs);
1423ca3ac6ddSLeila Ghaffari                                   CHKERRQ(ierr);
1424*1e150236SLeila Ghaffari   }
1425*1e150236SLeila Ghaffari   if (qf_ifunctionSur) {
1426*1e150236SLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, user->op_ifunction_vol, qf_ifunctionSur, qf_setupSur, height,
1427ca3ac6ddSLeila Ghaffari                                   bc.noutflow, bc.outflow, numP_Sur, numQ_Sur, qdatasizeSur,
1428*1e150236SLeila Ghaffari                                   NqptsSur, basisxSur, basisqSur, &user->op_ifunction);
1429ca3ac6ddSLeila Ghaffari                                   CHKERRQ(ierr);
1430*1e150236SLeila Ghaffari   }
1431cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1432ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1433ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1434ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1435ccaff030SJeremy L Thompson     .CtauS = CtauS,
1436ccaff030SJeremy L Thompson     .strong_form = strong_form,
1437ccaff030SJeremy L Thompson     .stabilization = stab,
1438ccaff030SJeremy L Thompson   };
1439ccaff030SJeremy L Thompson   switch (problemChoice) {
1440ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1441ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1442ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1443ccaff030SJeremy L Thompson           sizeof ctxNS);
14446a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
14456a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
14466a0edaf9SLeila Ghaffari           sizeof ctxNS);
1447ccaff030SJeremy L Thompson     break;
1448ccaff030SJeremy L Thompson   case NS_ADVECTION:
1449ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1450ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1451ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1452ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1453ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
14546a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
14556a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
14566a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
14576a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1458ccaff030SJeremy L Thompson   }
1459ccaff030SJeremy L Thompson 
1460ccaff030SJeremy L Thompson   // Set up PETSc context
1461ccaff030SJeremy L Thompson   // Set up units structure
1462ccaff030SJeremy L Thompson   units->meter = meter;
1463ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1464ccaff030SJeremy L Thompson   units->second = second;
1465ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1466ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1467ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1468ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1469ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1470ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1471ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1472ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1473ccaff030SJeremy L Thompson 
1474ccaff030SJeremy L Thompson   // Set up user structure
1475ccaff030SJeremy L Thompson   user->comm = comm;
1476ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1477ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1478ccaff030SJeremy L Thompson   user->units = units;
1479ccaff030SJeremy L Thompson   user->dm = dm;
1480ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1481ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1482ccaff030SJeremy L Thompson   user->ceed = ceed;
1483ccaff030SJeremy L Thompson 
14848b982baeSLeila Ghaffari   // Calculate qdata and ICs
1485ccaff030SJeremy L Thompson   // Set up state global and local vectors
1486ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1487ccaff030SJeremy L Thompson 
1488cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1489ccaff030SJeremy L Thompson 
1490ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1491ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
14928b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
149384d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1494ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1495ccaff030SJeremy L Thompson 
149684d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
149784d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1498ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1499ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1500ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1501ccaff030SJeremy L Thompson     // slower execution.
1502ccaff030SJeremy L Thompson     Vec Qbc;
1503ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1504ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1505ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1506ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1507ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1508ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1509ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
151084d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
151184d34d69SLeila Ghaffari     CHKERRQ(ierr);
1512ccaff030SJeremy L Thompson   }
1513ccaff030SJeremy L Thompson 
1514ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1515ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1516ccaff030SJeremy L Thompson   // Gather initial Q values
1517ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1518ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1519ccaff030SJeremy L Thompson     PetscViewer viewer;
1520ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1521ccaff030SJeremy L Thompson     // Read input
1522ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1523ccaff030SJeremy L Thompson                          user->outputfolder);
1524ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1525ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1526ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1527ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1528ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1529ccaff030SJeremy L Thompson   }
1530ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1531ccaff030SJeremy L Thompson 
1532ccaff030SJeremy L Thompson // Create and setup TS
1533ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1534ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1535ccaff030SJeremy L Thompson   if (implicit) {
1536ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1537ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1538ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1539ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1540ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1541ccaff030SJeremy L Thompson     }
1542ccaff030SJeremy L Thompson   } else {
1543ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1544ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1545ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1546ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1547ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1548ccaff030SJeremy L Thompson   }
1549ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1550ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1551ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
155284d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1553ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1554ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1555ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1556ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1557ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1558ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
155984d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1560ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1561ccaff030SJeremy L Thompson     }
1562ccaff030SJeremy L Thompson   } else { // continue from time of last output
1563ccaff030SJeremy L Thompson     PetscReal time;
1564ccaff030SJeremy L Thompson     PetscInt count;
1565ccaff030SJeremy L Thompson     PetscViewer viewer;
1566ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1567ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1568ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1569ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1570ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1571ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1572ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1573ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1574ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1575ccaff030SJeremy L Thompson   }
157684d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1577ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1578ccaff030SJeremy L Thompson   }
1579ccaff030SJeremy L Thompson 
1580ccaff030SJeremy L Thompson   // Solve
1581ccaff030SJeremy L Thompson   start = MPI_Wtime();
1582ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1583ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1584ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1585ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1586ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1587ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
158884d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1589ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
159084d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1591ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1592ccaff030SJeremy L Thompson   }
1593ccaff030SJeremy L Thompson 
1594ccaff030SJeremy L Thompson   // Get error
159584d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1596ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1597ccaff030SJeremy L Thompson     PetscReal norm;
1598ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1599ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1600ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1601ccaff030SJeremy L Thompson 
160284d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
160384d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1604ccaff030SJeremy L Thompson 
1605ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1607cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1608ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1609ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1610ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
161184d34d69SLeila Ghaffari     // Clean up vectors
161284d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
161384d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1614ccaff030SJeremy L Thompson   }
1615ccaff030SJeremy L Thompson 
1616ccaff030SJeremy L Thompson   // Output Statistics
1617ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
161884d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1619ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1620ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1621ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1622ccaff030SJeremy L Thompson   }
162384d34d69SLeila Ghaffari   // Output numerical values from command line
162484d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
162584d34d69SLeila Ghaffari 
162684d34d69SLeila Ghaffari   // compare reference solution values with current run
162784d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
162884d34d69SLeila Ghaffari     PetscViewer viewer;
162984d34d69SLeila Ghaffari     // Read reference file
163084d34d69SLeila Ghaffari     Vec Qref;
163184d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
163284d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
163384d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
163484d34d69SLeila Ghaffari     CHKERRQ(ierr);
163584d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
163684d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
163784d34d69SLeila Ghaffari 
163884d34d69SLeila Ghaffari     // Compute error with respect to reference solution
163984d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
164084d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
164184d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
164284d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
164384d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
164484d34d69SLeila Ghaffari     // Check error
164584d34d69SLeila Ghaffari     if (error > test->testtol) {
164684d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
164784d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
164884d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
164984d34d69SLeila Ghaffari     }
165084d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
165184d34d69SLeila Ghaffari   }
16529cf88b28Svaleriabarra 
1653ccaff030SJeremy L Thompson   // Clean up libCEED
16548b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1655ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1656ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1657ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1658ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
165984d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
166084d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
166184d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
166284d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
166384d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
166484d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1665ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1666ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1667ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1668ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1669ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1670ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
16716a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
16726a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
16736a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
16746a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
16756a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
16766a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
16776a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
16786a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
16796a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
1680ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1681ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1682ccaff030SJeremy L Thompson 
1683ccaff030SJeremy L Thompson   // Clean up PETSc
1684ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1685ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1686ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1687ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1688ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1689ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1690ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1691ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1692ccaff030SJeremy L Thompson   return PetscFinalize();
1693ccaff030SJeremy L Thompson }
1694ccaff030SJeremy L Thompson 
1695