xref: /libCEED/examples/fluids/navierstokes.c (revision 32b5ec5fd58411688994efaf910dd0fa1abe4e10)
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;
2221e150236SLeila 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);
303*32b5ec5fSJed Brown     bool flip = false;
304*32b5ec5fSJed Brown     if (height > 0) {
305*32b5ec5fSJed Brown       PetscInt numCells, numFaces, start = -1;
306*32b5ec5fSJed Brown       const PetscInt *orients, *faces, *cells;
307*32b5ec5fSJed Brown       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
308*32b5ec5fSJed Brown       ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr);
309*32b5ec5fSJed Brown       if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Expected one cell in support of exterior face, but got %D cells", numCells);
310*32b5ec5fSJed Brown       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
311*32b5ec5fSJed Brown       ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr);
312*32b5ec5fSJed Brown       for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;}
313*32b5ec5fSJed Brown       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %D in cone of its support", c);
314*32b5ec5fSJed Brown       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
315*32b5ec5fSJed Brown       if (orients[start] < 0) flip = true;
316*32b5ec5fSJed Brown     }
31784d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
31884d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
31984d34d69SLeila Ghaffari           c);
320ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
321ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
322*32b5ec5fSJed Brown       PetscInt ii = i;
323*32b5ec5fSJed Brown       if (flip) {
324*32b5ec5fSJed Brown 	if (P == nnodes) ii = nnodes - 1 - i;
325*32b5ec5fSJed Brown 	else if (P*P == nnodes) {
326*32b5ec5fSJed Brown 	  PetscInt row = i / P, col = i % P;
327*32b5ec5fSJed Brown 	  ii = row + col * P;
328*32b5ec5fSJed Brown 	} else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with %D nodes != P (%D) or P^2", nnodes, P);
329*32b5ec5fSJed Brown       }
330ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
331ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
332ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
333ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
334*32b5ec5fSJed Brown           if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j])
335*32b5ec5fSJed Brown               != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j)
336ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
337ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
338*32b5ec5fSJed Brown                      c, ii, f, j);
339ccaff030SJeremy L Thompson         }
340ccaff030SJeremy L Thompson       }
341ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
342*32b5ec5fSJed Brown       PetscInt loc = Involute(indices[ii*ncomp[0]]);
3436f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
344ccaff030SJeremy L Thompson     }
34584d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
34684d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
34784d34d69SLeila Ghaffari     CHKERRQ(ierr);
348ccaff030SJeremy L Thompson   }
3490c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3500c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3510c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
352ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3530c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3540c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3550c6c0b13SLeila Ghaffari 
356ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
357ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
358ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3596f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3606f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3616f55dfd5Svaleriabarra                             Erestrict);
362ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
364ccaff030SJeremy L Thompson }
365ccaff030SJeremy L Thompson 
366c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3671e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
3681e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
3691e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
3701e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
371c96c872fSLeila Ghaffari 
372c96c872fSLeila Ghaffari   DM dmcoord;
3731e150236SLeila Ghaffari   CeedInt dim, localNelem;
3741e150236SLeila Ghaffari   CeedInt Qdim;
375c96c872fSLeila Ghaffari   PetscErrorCode ierr;
376c96c872fSLeila Ghaffari 
377c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
3781e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
3791e150236SLeila Ghaffari   dim -= height;
3801e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
381c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
382c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
383c96c872fSLeila Ghaffari   CHKERRQ(ierr);
384c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
385c96c872fSLeila Ghaffari   CHKERRQ(ierr);
386c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
387c96c872fSLeila Ghaffari   CHKERRQ(ierr);
388c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
389c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
390c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
391c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
392c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
393c96c872fSLeila Ghaffari }
394c96c872fSLeila Ghaffari 
3951e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
3961e150236SLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, CeedOperator op_applyVol,
397ca3ac6ddSLeila Ghaffari     CeedQFunction qf_applySur, CeedQFunction qf_setupSur, CeedInt height, PetscInt nFace,
398ca3ac6ddSLeila Ghaffari     PetscInt value[6], CeedInt numP_Sur, CeedInt numQ_Sur, CeedInt qdatasizeSur, CeedInt NqptsSur,
3991e150236SLeila Ghaffari     CeedBasis basisxSur, CeedBasis basisqSur, CeedOperator *op_apply) {
400ca3ac6ddSLeila Ghaffari 
401ca3ac6ddSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
402ca3ac6ddSLeila Ghaffari   PetscInt dim, localNelemSur[6];
4031e150236SLeila Ghaffari   Vec Xloc;
4041e150236SLeila Ghaffari   CeedVector xcorners, qdataSur[6];
405ca3ac6ddSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
406ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
4071e150236SLeila Ghaffari   PetscScalar *x;
4081e150236SLeila Ghaffari   PetscInt lsize;
409ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
410ca3ac6ddSLeila Ghaffari 
411ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
412ca3ac6ddSLeila Ghaffari   // Composite Operaters
413ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
4141e150236SLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_applyVol); // Apply a Sub-Operator for the volume
415ca3ac6ddSLeila Ghaffari 
416ca3ac6ddSLeila Ghaffari   if (nFace) {
4171e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
4181e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
4191e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
4201e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
4211e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
4221e150236SLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
423ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
4241e150236SLeila Ghaffari     // Create CEED Operator for each In/OutFlow faces
425ca3ac6ddSLeila Ghaffari     for(CeedInt i=0; i<nFace; i++) {
4261e150236SLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, value[i], numP_Sur,
427ca3ac6ddSLeila Ghaffari                                      numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
428ca3ac6ddSLeila Ghaffari                                      &restrictqdiSur[i]); CHKERRQ(ierr);
4291e150236SLeila Ghaffari 
4301e150236SLeila Ghaffari       // Create the CEED vectors that will be needed in boundary setup
431ca3ac6ddSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
432ca3ac6ddSLeila Ghaffari       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
433ca3ac6ddSLeila Ghaffari 
4341e150236SLeila Ghaffari       // Create the operator that builds the quadrature data for the In/OutFlow operator
435ca3ac6ddSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
436ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, CEED_VECTOR_ACTIVE);
437ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
438ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
439ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
440ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4411e150236SLeila Ghaffari 
4421e150236SLeila Ghaffari       // Create In/OutFlow operator
443ca3ac6ddSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
444ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
445ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
446ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
447ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, xcorners);
448ca3ac6ddSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
449ca3ac6ddSLeila Ghaffari 
4501e150236SLeila Ghaffari       // Apply CEED operator for boundary setup
4511e150236SLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i], CEED_REQUEST_IMMEDIATE);
4521e150236SLeila Ghaffari 
4531e150236SLeila Ghaffari       // Apply Sub-Operator for In/OutFlow BCs
454ca3ac6ddSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
455ca3ac6ddSLeila Ghaffari     }
4561e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
457ca3ac6ddSLeila Ghaffari   }
458ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
459ca3ac6ddSLeila Ghaffari }
460ca3ac6ddSLeila Ghaffari 
461ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
462ccaff030SJeremy L Thompson   PetscErrorCode ierr;
463ccaff030SJeremy L Thompson   PetscInt m;
464ccaff030SJeremy L Thompson 
465ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
466ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
468ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
469ccaff030SJeremy L Thompson }
470ccaff030SJeremy L Thompson 
471ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
472ccaff030SJeremy L Thompson   PetscErrorCode ierr;
473ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
474ccaff030SJeremy L Thompson   PetscScalar *a;
475ccaff030SJeremy L Thompson 
476ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
477ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
479ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
48084d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
48184d34d69SLeila Ghaffari                                   mpetsc, mceed);
482ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
483ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
484ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
485ccaff030SJeremy L Thompson }
486ccaff030SJeremy L Thompson 
487ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
488ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
489ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
490ccaff030SJeremy L Thompson   PetscErrorCode ierr;
491ccaff030SJeremy L Thompson   Vec Qbc;
492ccaff030SJeremy L Thompson 
493ccaff030SJeremy L Thompson   PetscFunctionBegin;
494ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
495ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
496ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
498ccaff030SJeremy L Thompson }
499ccaff030SJeremy L Thompson 
500ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
501ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
502ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
503ccaff030SJeremy L Thompson   PetscErrorCode ierr;
504ccaff030SJeremy L Thompson   User user = *(User *)userData;
505ccaff030SJeremy L Thompson   PetscScalar *q, *g;
506ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
507ccaff030SJeremy L Thompson 
508ccaff030SJeremy L Thompson   // Global-to-local
509ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
510ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
511ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
512ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
513ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
514ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
515ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
516ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
517ccaff030SJeremy L Thompson 
518ccaff030SJeremy L Thompson   // Ceed Vectors
519ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
520ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
521ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
522ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
523ccaff030SJeremy L Thompson 
524ccaff030SJeremy L Thompson   // Apply CEED operator
525ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
526ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
527ccaff030SJeremy L Thompson 
528ccaff030SJeremy L Thompson   // Restore vectors
529ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
530ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
531ccaff030SJeremy L Thompson 
532ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson 
535ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
536ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
537ccaff030SJeremy L Thompson   CHKERRQ(ierr);
538ccaff030SJeremy L Thompson 
539ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
542ccaff030SJeremy L Thompson }
543ccaff030SJeremy L Thompson 
544ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
545ccaff030SJeremy L Thompson                                    void *userData) {
546ccaff030SJeremy L Thompson   PetscErrorCode ierr;
547ccaff030SJeremy L Thompson   User user = *(User *)userData;
548ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
549ccaff030SJeremy L Thompson   PetscScalar *g;
550ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
551ccaff030SJeremy L Thompson 
552ccaff030SJeremy L Thompson   // Global-to-local
553ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
554ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
555ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
556ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
557ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
560ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
562ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
563ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson 
565ccaff030SJeremy L Thompson   // Ceed Vectors
566ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
567ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
568ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
569ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
570ccaff030SJeremy L Thompson                      (PetscScalar *)q);
571ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
572ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
573ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
574ccaff030SJeremy L Thompson 
575ccaff030SJeremy L Thompson   // Apply CEED operator
576ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
577ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
578ccaff030SJeremy L Thompson 
579ccaff030SJeremy L Thompson   // Restore vectors
580ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
581ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
582ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson 
584ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
585ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
586ccaff030SJeremy L Thompson 
587ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
589ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
591ccaff030SJeremy L Thompson }
592ccaff030SJeremy L Thompson 
593ccaff030SJeremy L Thompson // User provided TS Monitor
594ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
595ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
596ccaff030SJeremy L Thompson   User user = ctx;
597ccaff030SJeremy L Thompson   Vec Qloc;
598ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
599ccaff030SJeremy L Thompson   PetscViewer viewer;
600ccaff030SJeremy L Thompson   PetscErrorCode ierr;
601ccaff030SJeremy L Thompson 
602ccaff030SJeremy L Thompson   // Set up output
603ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
604ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
605ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
606ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
607ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
609ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
610ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson   // Output
613ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
614ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
615ccaff030SJeremy L Thompson   CHKERRQ(ierr);
616ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
617ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
618ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
6199d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
620ccaff030SJeremy L Thompson   if (user->dmviz) {
621ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
622ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
623ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
624ccaff030SJeremy L Thompson 
625ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
626ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
627ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
628ccaff030SJeremy L Thompson     CHKERRQ(ierr);
629ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
630ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
631ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
632ccaff030SJeremy L Thompson     CHKERRQ(ierr);
633ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
634ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
635ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
636ccaff030SJeremy L Thompson     CHKERRQ(ierr);
637ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
638ccaff030SJeremy L Thompson                               filepath_refined,
639ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
640ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
641ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
643ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
644ccaff030SJeremy L Thompson   }
645ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
646ccaff030SJeremy L Thompson 
647ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
648ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
649ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
650ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
651ccaff030SJeremy L Thompson   CHKERRQ(ierr);
652ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
653ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
654ccaff030SJeremy L Thompson 
655ccaff030SJeremy L Thompson   // Save time stamp
656ccaff030SJeremy L Thompson   // Dimensionalize time back
657ccaff030SJeremy L Thompson   time /= user->units->second;
658ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
659ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
660ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
661ccaff030SJeremy L Thompson   CHKERRQ(ierr);
662ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
663ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
664ccaff030SJeremy L Thompson   #else
665ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
666ccaff030SJeremy L Thompson   #endif
667ccaff030SJeremy L Thompson   CHKERRQ(ierr);
668ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
669ccaff030SJeremy L Thompson 
670ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
671ccaff030SJeremy L Thompson }
672ccaff030SJeremy L Thompson 
67384d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
674ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
675ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
676ccaff030SJeremy L Thompson   PetscErrorCode ierr;
677ccaff030SJeremy L Thompson   CeedVector multlvec;
678ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
679ccaff030SJeremy L Thompson 
680ccaff030SJeremy L Thompson   ctxSetup->time = time;
681ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
682ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
683ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
684ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
685ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
686ccaff030SJeremy L Thompson 
687ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
688ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
689ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
690ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
691ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
692ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
693ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
694ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
695ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
696ccaff030SJeremy L Thompson   CHKERRQ(ierr);
697ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
698ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
699ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
700ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
701ccaff030SJeremy L Thompson 
702ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
703ccaff030SJeremy L Thompson }
704ccaff030SJeremy L Thompson 
705ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
706ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
707ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
708ccaff030SJeremy L Thompson   PetscErrorCode ierr;
709ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
710ccaff030SJeremy L Thompson   CeedOperator op_mass;
711ccaff030SJeremy L Thompson   CeedVector mceed;
712ccaff030SJeremy L Thompson   Vec Mloc;
713ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
714ccaff030SJeremy L Thompson 
715ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
716ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
717ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
718ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
719ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
720ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
721ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
722ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
723ccaff030SJeremy L Thompson 
724ccaff030SJeremy L Thompson   // Create the mass operator
725ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
726ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
727ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
728ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
729ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
730ccaff030SJeremy L Thompson 
731ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
732ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
733ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
734ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
735ccaff030SJeremy L Thompson 
736ccaff030SJeremy L Thompson   {
737ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
738ccaff030SJeremy L Thompson     CeedVector onesvec;
739ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
740ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
741ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
742ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
743ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
744ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
745ccaff030SJeremy L Thompson   }
746ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
747ccaff030SJeremy L Thompson 
748ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
749ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
750ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
751ccaff030SJeremy L Thompson 
752ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
753ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
754ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
755ccaff030SJeremy L Thompson }
756ccaff030SJeremy L Thompson 
75784d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
758ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
759ccaff030SJeremy L Thompson   PetscErrorCode ierr;
760ccaff030SJeremy L Thompson 
761ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
762ccaff030SJeremy L Thompson   {
763ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
764ccaff030SJeremy L Thompson     PetscFE fe;
765ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
766ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
767ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
76832ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
769ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
770ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
771ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
77207af6069Svaleriabarra     {
77307af6069Svaleriabarra       PetscInt comps[1] = {1};
77407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
77507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
77607af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
77707af6069Svaleriabarra       comps[0] = 2;
77807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
77907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
78007af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
78107af6069Svaleriabarra       comps[0] = 3;
78207af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
78307af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
78407af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
78507af6069Svaleriabarra     }
78684d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
78784d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
78884d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
78984d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
79084d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
79184d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
79284d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
79384d34d69SLeila Ghaffari 
79484d34d69SLeila Ghaffari           }
79584d34d69SLeila Ghaffari         }
79684d34d69SLeila Ghaffari       }
79784d34d69SLeila Ghaffari     }
79884d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
79984d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
80084d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
80184d34d69SLeila Ghaffari     {
80284d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
80384d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
80484d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
80584d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
80684d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
80784d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
80884d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
80984d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
81084d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
81184d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
81284d34d69SLeila Ghaffari       } else
81384d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
81484d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
81584d34d69SLeila Ghaffari     }
816ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
817ccaff030SJeremy L Thompson     CHKERRQ(ierr);
818ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
819ccaff030SJeremy L Thompson   }
820ccaff030SJeremy L Thompson   {
821ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
822ccaff030SJeremy L Thompson     PetscSection section;
823ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
824ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
825ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
826ccaff030SJeremy L Thompson     CHKERRQ(ierr);
827ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
828ccaff030SJeremy L Thompson     CHKERRQ(ierr);
829ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
830ccaff030SJeremy L Thompson     CHKERRQ(ierr);
831ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
832ccaff030SJeremy L Thompson     CHKERRQ(ierr);
833ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
834ccaff030SJeremy L Thompson     CHKERRQ(ierr);
835ccaff030SJeremy L Thompson   }
836ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
837ccaff030SJeremy L Thompson }
838ccaff030SJeremy L Thompson 
839ccaff030SJeremy L Thompson int main(int argc, char **argv) {
840ccaff030SJeremy L Thompson   PetscInt ierr;
841ccaff030SJeremy L Thompson   MPI_Comm comm;
84284d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
843ccaff030SJeremy L Thompson   Mat interpviz;
844ccaff030SJeremy L Thompson   TS ts;
845ccaff030SJeremy L Thompson   TSAdapt adapt;
846ccaff030SJeremy L Thompson   User user;
847ccaff030SJeremy L Thompson   Units units;
848ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
84984d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
850ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
851ccaff030SJeremy L Thompson   PetscMPIInt rank;
852ccaff030SJeremy L Thompson   PetscScalar ftime;
853ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
854ccaff030SJeremy L Thompson   Ceed ceed;
85584d34d69SLeila Ghaffari   CeedInt numP, numQ;
856cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
85784d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
85884d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
859cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
860cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
861ccaff030SJeremy L Thompson   CeedScalar Rd;
86284d34d69SLeila Ghaffari   CeedMemType memtyperequested;
863ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
864ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
865ccaff030SJeremy L Thompson   problemType problemChoice;
866ccaff030SJeremy L Thompson   problemData *problem = NULL;
867ccaff030SJeremy L Thompson   StabilizationType stab;
86884d34d69SLeila Ghaffari   testType testChoice;
86984d34d69SLeila Ghaffari   testData *test = NULL;
87084d34d69SLeila Ghaffari   PetscBool implicit;
871cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
872ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
87384d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
87484d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
875ccaff030SJeremy L Thompson   };
876ccaff030SJeremy L Thompson   double start, cpu_time_used;
87784d34d69SLeila Ghaffari   // Check PETSc CUDA support
87884d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
87984d34d69SLeila Ghaffari   // *INDENT-OFF*
88084d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
88184d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
88284d34d69SLeila Ghaffari   #else
88384d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
88484d34d69SLeila Ghaffari   #endif
88584d34d69SLeila Ghaffari   // *INDENT-ON*
886ccaff030SJeremy L Thompson 
887ccaff030SJeremy L Thompson   // Create the libCEED contexts
888ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
889ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
890ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
891ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
892ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
893ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
894ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
895ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
896ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
897ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
898ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
899ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
900ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
901ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
902ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
903ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
904ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
905ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
906ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
907ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
908ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
909ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
910ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
911ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
912ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
913ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
91484d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
91584d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
916ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
917ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
918ccaff030SJeremy L Thompson 
919ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
920ccaff030SJeremy L Thompson   if (ierr) return ierr;
921ccaff030SJeremy L Thompson 
922ccaff030SJeremy L Thompson   // Allocate PETSc context
923ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
924ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
925ccaff030SJeremy L Thompson 
926ccaff030SJeremy L Thompson   // Parse command line options
927ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
928ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
929ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
930ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
931ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
932ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
93384d34d69SLeila Ghaffari   testChoice = TEST_NONE;
93484d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
93584d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
93684d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
93784d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
93884d34d69SLeila Ghaffari   test = &testOptions[testChoice];
939ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
940ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
941ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
942ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
943ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
944ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
945ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
946ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
947ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
948ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
949ccaff030SJeremy L Thompson   CHKERRQ(ierr);
95084d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
95184d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
95284d34d69SLeila Ghaffari     CHKERRQ(ierr);
95384d34d69SLeila Ghaffari   }
954ccaff030SJeremy L Thompson   {
9557573aee6SJed Brown     PetscInt len;
9567573aee6SJed Brown     PetscBool flg;
9577573aee6SJed Brown     ierr = PetscOptionsIntArray("-bc_outflow",
9587573aee6SJed Brown                               "Use outflow boundary conditions on this list of faces",
9597573aee6SJed Brown                               NULL, bc.outflow,
9607573aee6SJed Brown                               (len = sizeof(bc.outflow) / sizeof(bc.outflow[0]),
9617573aee6SJed Brown                               &len), &flg); CHKERRQ(ierr);
9627573aee6SJed Brown     if (flg) {
9637573aee6SJed Brown       bc.noutflow = len;
9647573aee6SJed Brown       // Using outflow boundaries disables automatic wall/slip boundaries (they must be set explicitly)
9657573aee6SJed Brown       bc.nwall = 0;
9667573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
9677573aee6SJed Brown     }
968ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
969ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
970ccaff030SJeremy L Thompson                                 NULL, bc.walls,
971ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
972ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
9737573aee6SJed Brown     if (flg) {
9747573aee6SJed Brown       bc.nwall = len;
9757573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
9767573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
9777573aee6SJed Brown     }
978ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
979ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
980ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
981ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
982ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
983ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
984ccaff030SJeremy L Thompson                                    &len), &flg);
985ccaff030SJeremy L Thompson       CHKERRQ(ierr);
98684d34d69SLeila Ghaffari       if (flg) {
98784d34d69SLeila Ghaffari         bc.nslip[j] = len;
98884d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
98984d34d69SLeila Ghaffari       }
990ccaff030SJeremy L Thompson     }
991ccaff030SJeremy L Thompson   }
992cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
993cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
994cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
995ccaff030SJeremy L Thompson   CHKERRQ(ierr);
996ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
997ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
998ccaff030SJeremy L Thompson   meter = fabs(meter);
999ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1000ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   second = fabs(second);
1002ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1003ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1004ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
1005ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
1006ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
1007ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1008ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
1009ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
1010ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
1011ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
1012ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
1013ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
1014ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
1015ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
1016ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
1017ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1018ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1019ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1020ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1021ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1022ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1023ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1024ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1025ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1026ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1027ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1028ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1029ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1030ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1031ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1032ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
103384d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
103484d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
103584d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
103684d34d69SLeila Ghaffari     CHKERRQ(ierr);
103784d34d69SLeila Ghaffari   }
1038ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1039ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1040ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1041ccaff030SJeremy L Thompson   CHKERRQ(ierr);
104284d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
104384d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
104484d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
104584d34d69SLeila Ghaffari     CHKERRQ(ierr);
104684d34d69SLeila Ghaffari   }
1047ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1048ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1049ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1050ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1051ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1052ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1053ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1054ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1055ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1056ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1057ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1058ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1059ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1060ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
1061ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
1062ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1063ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1064ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1065ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1066ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1067ccaff030SJeremy L Thompson   n = problem->dim;
1068ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1069ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1070ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1071ccaff030SJeremy L Thompson   {
1072ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1073ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1074ccaff030SJeremy L Thompson     if (norm > 0) {
1075ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1076ccaff030SJeremy L Thompson     }
1077ccaff030SJeremy L Thompson   }
1078ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1079ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1080ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1081ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1082ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
108384d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
108484d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
108584d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
108684d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
108784d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1088ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1089ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1090ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
109184d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
109284d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
109384d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
109484d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
109584d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
109684d34d69SLeila Ghaffari   CHKERRQ(ierr);
1097ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1098ccaff030SJeremy L Thompson 
1099ccaff030SJeremy L Thompson   // Define derived units
1100ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1101ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1102ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1103ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1104ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1105ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1106ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
1107ccaff030SJeremy L Thompson 
1108ccaff030SJeremy L Thompson   // Scale variables to desired units
1109ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1110ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1111ccaff030SJeremy L Thompson   P0 *= Pascal;
1112ccaff030SJeremy L Thompson   N *= (1./second);
1113ccaff030SJeremy L Thompson   cv *= JperkgK;
1114ccaff030SJeremy L Thompson   cp *= JperkgK;
1115ccaff030SJeremy L Thompson   Rd = cp - cv;
1116ccaff030SJeremy L Thompson   g *= mpersquareds;
1117ccaff030SJeremy L Thompson   mu *= Pascal * second;
1118ccaff030SJeremy L Thompson   k *= WpermK;
1119ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1120ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1121ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1122ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1123ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1124ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1125ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1126ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1127ccaff030SJeremy L Thompson 
1128ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1129cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1130ccaff030SJeremy L Thompson   // Set up the libCEED context
1131ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1132ccaff030SJeremy L Thompson     .theta0 = theta0,
1133ccaff030SJeremy L Thompson     .thetaC = thetaC,
1134ccaff030SJeremy L Thompson     .P0 = P0,
1135ccaff030SJeremy L Thompson     .N = N,
1136ccaff030SJeremy L Thompson     .cv = cv,
1137ccaff030SJeremy L Thompson     .cp = cp,
1138ccaff030SJeremy L Thompson     .Rd = Rd,
1139ccaff030SJeremy L Thompson     .g = g,
1140ccaff030SJeremy L Thompson     .rc = rc,
1141ccaff030SJeremy L Thompson     .lx = lx,
1142ccaff030SJeremy L Thompson     .ly = ly,
1143ccaff030SJeremy L Thompson     .lz = lz,
1144ccaff030SJeremy L Thompson     .center[0] = center[0],
1145ccaff030SJeremy L Thompson     .center[1] = center[1],
1146ccaff030SJeremy L Thompson     .center[2] = center[2],
1147ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1148ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1149ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
1150ccaff030SJeremy L Thompson     .time = 0,
1151ccaff030SJeremy L Thompson   };
1152ccaff030SJeremy L Thompson 
115384d34d69SLeila Ghaffari   // Create the mesh
1154ccaff030SJeremy L Thompson   {
1155ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1156ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
115784d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1158ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1159ccaff030SJeremy L Thompson   }
116084d34d69SLeila Ghaffari 
116184d34d69SLeila Ghaffari   // Distribute the mesh over processes
116284d34d69SLeila Ghaffari   {
1163ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1164ccaff030SJeremy L Thompson     PetscPartitioner part;
1165ccaff030SJeremy L Thompson 
1166ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1167ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1168ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1169ccaff030SJeremy L Thompson     if (dmDist) {
1170ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1171ccaff030SJeremy L Thompson       dm  = dmDist;
1172ccaff030SJeremy L Thompson     }
1173ccaff030SJeremy L Thompson   }
1174ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1175ccaff030SJeremy L Thompson 
117684d34d69SLeila Ghaffari   // Setup DM
1177ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1178ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
117984d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
118084d34d69SLeila Ghaffari 
118184d34d69SLeila Ghaffari   // Refine DM for high-order viz
1182ccaff030SJeremy L Thompson   dmviz = NULL;
1183ccaff030SJeremy L Thompson   interpviz = NULL;
1184ccaff030SJeremy L Thompson   if (viz_refine) {
1185ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1186ff6701fcSJed Brown 
1187ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1188ff6701fcSJed Brown     dmhierarchy[0] = dm;
118984d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1190ff6701fcSJed Brown       Mat interp_next;
1191ff6701fcSJed Brown 
1192ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1193ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1194ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1195ff6701fcSJed Brown       d = (d + 1) / 2;
1196ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1197ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1198ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1199ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1200ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1201ff6701fcSJed Brown       else {
1202ff6701fcSJed Brown         Mat C;
1203ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1204ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1205ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1206ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1207ff6701fcSJed Brown         interpviz = C;
1208ff6701fcSJed Brown       }
1209ff6701fcSJed Brown     }
1210cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1211ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1212cb3e2689Svaleriabarra     }
1213ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1214ccaff030SJeremy L Thompson   }
1215ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1216ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1217ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1218ccaff030SJeremy L Thompson   lnodes /= ncompq;
1219ccaff030SJeremy L Thompson 
122084d34d69SLeila Ghaffari   // Initialize CEED
122184d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
122284d34d69SLeila Ghaffari   // Set memtype
122384d34d69SLeila Ghaffari   CeedMemType memtypebackend;
122484d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
122584d34d69SLeila Ghaffari   // Check memtype compatibility
122684d34d69SLeila Ghaffari   if (!setmemtyperequest)
122784d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
122884d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
122984d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
123084d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
123184d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
123284d34d69SLeila Ghaffari 
123384d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
123484d34d69SLeila Ghaffari   numP = degree + 1;
123584d34d69SLeila Ghaffari   numQ = numP + qextra;
123684d34d69SLeila Ghaffari 
123784d34d69SLeila Ghaffari     // Print summary
123884d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1239ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1240ccaff030SJeremy L Thompson     int comm_size;
1241ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1242ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1243ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
124484d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1245ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1246ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1247ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
124884d34d69SLeila Ghaffari     const char *usedresource;
124984d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1250ccaff030SJeremy L Thompson 
125184d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
125284d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
125384d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
125484d34d69SLeila Ghaffari                        "  Problem:\n"
125584d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
125684d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
125784d34d69SLeila Ghaffari                        "  PETSc:\n"
125884d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
125984d34d69SLeila Ghaffari                        "  libCEED:\n"
126084d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
126184d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
126284d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
126384d34d69SLeila Ghaffari                        "  Mesh:\n"
126484d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
126584d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
126684d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
126784d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
126884d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
126984d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
127084d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
127184d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
127284d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
127384d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
127484d34d69SLeila Ghaffari                        (setmemtyperequest) ?
127584d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
127684d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
127784d34d69SLeila Ghaffari     CHKERRQ(ierr);
12780c6c0b13SLeila Ghaffari   }
12790c6c0b13SLeila Ghaffari 
1280ccaff030SJeremy L Thompson   // Set up global mass vector
1281ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1282ccaff030SJeremy L Thompson 
128384d34d69SLeila Ghaffari   // Set up libCEED
1284ccaff030SJeremy L Thompson   // CEED Bases
1285ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
128684d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
128784d34d69SLeila Ghaffari                                   &basisq);
128884d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
128984d34d69SLeila Ghaffari                                   &basisx);
129084d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
129184d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1292ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1293ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1294ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1295ccaff030SJeremy L Thompson 
1296ccaff030SJeremy L Thompson   // CEED Restrictions
12971e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
129884d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
129984d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1300ccaff030SJeremy L Thompson 
1301ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1302ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1303ccaff030SJeremy L Thompson 
1304ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1305bd910870SLeila Ghaffari   CeedInt NqptsVol;
130684d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
130784d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
13088b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
130984d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1310ccaff030SJeremy L Thompson 
1311ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1312ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1313ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1314ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1315ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
13168b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1317ccaff030SJeremy L Thompson 
1318ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1319ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1320ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1321ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1322ccaff030SJeremy L Thompson 
1323ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1324ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1325ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1326ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1327ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1328ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13298b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1330ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1331ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1332ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1333ccaff030SJeremy L Thompson   }
1334ccaff030SJeremy L Thompson 
1335ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1336ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1337ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1338ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1339ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1340ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1341ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
13428b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1343ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1344ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1345ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1346ccaff030SJeremy L Thompson   }
1347ccaff030SJeremy L Thompson 
1348ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1349ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
135084d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1351ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
135284d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
135384d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1354ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1355ccaff030SJeremy L Thompson 
1356ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1357ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
135884d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
135984d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1360ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1361ccaff030SJeremy L Thompson 
136284d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
136384d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
136484d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1365ccaff030SJeremy L Thompson 
1366ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1367ccaff030SJeremy L Thompson     CeedOperator op;
1368ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
136984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
137084d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
137184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
13728b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
137384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
137484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
137584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1376d3630711SJed Brown     user->op_rhs_vol = op;
1377ccaff030SJeremy L Thompson   }
1378ccaff030SJeremy L Thompson 
1379ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1380ccaff030SJeremy L Thompson     CeedOperator op;
1381ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
138284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
138384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
138484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
138584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
13868b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
138784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
138884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
138984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1390d3630711SJed Brown     user->op_ifunction_vol = op;
1391ccaff030SJeremy L Thompson   }
1392ccaff030SJeremy L Thompson 
1393cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1394cfa64770SLeila Ghaffari   // Outflow Boundary Condition
13956a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
13966a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
13976a0edaf9SLeila Ghaffari   CeedInt height = 1;
13986a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
13991e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
14001e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1401cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1402cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1403cfa64770SLeila Ghaffari   CeedInt NqptsSur;
1404cfa64770SLeila Ghaffari   CeedQFunction qf_setupSur, qf_rhsSur, qf_ifunctionSur;
1405cfa64770SLeila Ghaffari 
1406cfa64770SLeila Ghaffari   // CEED bases for the boundaries
14076a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
14086a0edaf9SLeila Ghaffari                                   &basisqSur);
14096a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
14106a0edaf9SLeila Ghaffari                                   &basisxSur);
14116a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
14126a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
14136a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
14146a0edaf9SLeila Ghaffari 
1415cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
14166a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
14176a0edaf9SLeila Ghaffari                               &qf_setupSur);
14186a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
14196a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
14206a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14216a0edaf9SLeila Ghaffari 
14226a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
1423cfa64770SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator on the Surface
14246a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
14256a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
14266a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
14276a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14286a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
14296a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
14306a0edaf9SLeila Ghaffari   }
14316a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
14326a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
14336a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
14346a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
14356a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
14366a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14376a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
14386a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
14396a0edaf9SLeila Ghaffari   }
1440cfa64770SLeila Ghaffari   // Create CEED Operator for each face
14411e150236SLeila Ghaffari   if (qf_rhsSur) {
14421e150236SLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, user->op_rhs_vol, qf_rhsSur, qf_setupSur, height,
1443ca3ac6ddSLeila Ghaffari                                   bc.noutflow, bc.outflow, numP_Sur, numQ_Sur, qdatasizeSur,
14441e150236SLeila Ghaffari                                   NqptsSur, basisxSur, basisqSur, &user->op_rhs);
1445ca3ac6ddSLeila Ghaffari                                   CHKERRQ(ierr);
14461e150236SLeila Ghaffari   }
14471e150236SLeila Ghaffari   if (qf_ifunctionSur) {
14481e150236SLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, user->op_ifunction_vol, qf_ifunctionSur, qf_setupSur, height,
1449ca3ac6ddSLeila Ghaffari                                   bc.noutflow, bc.outflow, numP_Sur, numQ_Sur, qdatasizeSur,
14501e150236SLeila Ghaffari                                   NqptsSur, basisxSur, basisqSur, &user->op_ifunction);
1451ca3ac6ddSLeila Ghaffari                                   CHKERRQ(ierr);
14521e150236SLeila Ghaffari   }
1453cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1454ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1455ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1456ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1457ccaff030SJeremy L Thompson     .CtauS = CtauS,
1458ccaff030SJeremy L Thompson     .strong_form = strong_form,
1459ccaff030SJeremy L Thompson     .stabilization = stab,
1460ccaff030SJeremy L Thompson   };
1461ccaff030SJeremy L Thompson   switch (problemChoice) {
1462ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1463ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1464ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1465ccaff030SJeremy L Thompson           sizeof ctxNS);
14666a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
14676a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
14686a0edaf9SLeila Ghaffari           sizeof ctxNS);
1469ccaff030SJeremy L Thompson     break;
1470ccaff030SJeremy L Thompson   case NS_ADVECTION:
1471ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1472ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1473ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1474ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1475ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
14766a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
14776a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
14786a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
14796a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1480ccaff030SJeremy L Thompson   }
1481ccaff030SJeremy L Thompson 
1482ccaff030SJeremy L Thompson   // Set up PETSc context
1483ccaff030SJeremy L Thompson   // Set up units structure
1484ccaff030SJeremy L Thompson   units->meter = meter;
1485ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1486ccaff030SJeremy L Thompson   units->second = second;
1487ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1488ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1489ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1490ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1491ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1492ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1493ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1494ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1495ccaff030SJeremy L Thompson 
1496ccaff030SJeremy L Thompson   // Set up user structure
1497ccaff030SJeremy L Thompson   user->comm = comm;
1498ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1499ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1500ccaff030SJeremy L Thompson   user->units = units;
1501ccaff030SJeremy L Thompson   user->dm = dm;
1502ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1503ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1504ccaff030SJeremy L Thompson   user->ceed = ceed;
1505ccaff030SJeremy L Thompson 
15068b982baeSLeila Ghaffari   // Calculate qdata and ICs
1507ccaff030SJeremy L Thompson   // Set up state global and local vectors
1508ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1509ccaff030SJeremy L Thompson 
1510cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1511ccaff030SJeremy L Thompson 
1512ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1513ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
15148b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
151584d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1516ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1517ccaff030SJeremy L Thompson 
151884d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
151984d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1520ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1521ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1522ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1523ccaff030SJeremy L Thompson     // slower execution.
1524ccaff030SJeremy L Thompson     Vec Qbc;
1525ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1526ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1527ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1528ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1529ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1530ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1531ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
153284d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
153384d34d69SLeila Ghaffari     CHKERRQ(ierr);
1534ccaff030SJeremy L Thompson   }
1535ccaff030SJeremy L Thompson 
1536ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1537ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1538ccaff030SJeremy L Thompson   // Gather initial Q values
1539ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1540ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1541ccaff030SJeremy L Thompson     PetscViewer viewer;
1542ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1543ccaff030SJeremy L Thompson     // Read input
1544ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1545ccaff030SJeremy L Thompson                          user->outputfolder);
1546ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1547ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1548ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1549ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1550ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1551ccaff030SJeremy L Thompson   }
1552ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1553ccaff030SJeremy L Thompson 
1554ccaff030SJeremy L Thompson // Create and setup TS
1555ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1556ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1557ccaff030SJeremy L Thompson   if (implicit) {
1558ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1559ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1560ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1561ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1562ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1563ccaff030SJeremy L Thompson     }
1564ccaff030SJeremy L Thompson   } else {
1565ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1566ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1567ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1568ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1569ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1570ccaff030SJeremy L Thompson   }
1571ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1572ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1573ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
157484d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1575ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1576ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1577ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1578ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1579ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1580ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
158184d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1582ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1583ccaff030SJeremy L Thompson     }
1584ccaff030SJeremy L Thompson   } else { // continue from time of last output
1585ccaff030SJeremy L Thompson     PetscReal time;
1586ccaff030SJeremy L Thompson     PetscInt count;
1587ccaff030SJeremy L Thompson     PetscViewer viewer;
1588ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1589ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1590ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1591ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1592ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1593ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1594ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1595ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1596ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1597ccaff030SJeremy L Thompson   }
159884d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1599ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1600ccaff030SJeremy L Thompson   }
1601ccaff030SJeremy L Thompson 
1602ccaff030SJeremy L Thompson   // Solve
1603ccaff030SJeremy L Thompson   start = MPI_Wtime();
1604ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1605ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1607ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1608ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1609ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
161084d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1611ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
161284d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1613ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1614ccaff030SJeremy L Thompson   }
1615ccaff030SJeremy L Thompson 
1616ccaff030SJeremy L Thompson   // Get error
161784d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1618ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1619ccaff030SJeremy L Thompson     PetscReal norm;
1620ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1621ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1622ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1623ccaff030SJeremy L Thompson 
162484d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
162584d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1626ccaff030SJeremy L Thompson 
1627ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1628ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1629cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1630ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1631ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1632ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
163384d34d69SLeila Ghaffari     // Clean up vectors
163484d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
163584d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1636ccaff030SJeremy L Thompson   }
1637ccaff030SJeremy L Thompson 
1638ccaff030SJeremy L Thompson   // Output Statistics
1639ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
164084d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1641ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1642ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1643ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1644ccaff030SJeremy L Thompson   }
164584d34d69SLeila Ghaffari   // Output numerical values from command line
164684d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
164784d34d69SLeila Ghaffari 
164884d34d69SLeila Ghaffari   // compare reference solution values with current run
164984d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
165084d34d69SLeila Ghaffari     PetscViewer viewer;
165184d34d69SLeila Ghaffari     // Read reference file
165284d34d69SLeila Ghaffari     Vec Qref;
165384d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
165484d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
165584d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
165684d34d69SLeila Ghaffari     CHKERRQ(ierr);
165784d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
165884d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
165984d34d69SLeila Ghaffari 
166084d34d69SLeila Ghaffari     // Compute error with respect to reference solution
166184d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
166284d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
166384d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
166484d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
166584d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
166684d34d69SLeila Ghaffari     // Check error
166784d34d69SLeila Ghaffari     if (error > test->testtol) {
166884d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
166984d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
167084d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
167184d34d69SLeila Ghaffari     }
167284d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
167384d34d69SLeila Ghaffari   }
16749cf88b28Svaleriabarra 
1675ccaff030SJeremy L Thompson   // Clean up libCEED
16768b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1677ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1678ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1679ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1680ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
168184d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
168284d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
168384d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
168484d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
168584d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
168684d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1687ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1688ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1689ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1690ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1691ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1692ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
16936a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
16946a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
16956a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
16966a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
16976a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
16986a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
16996a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
17006a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
17016a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
1702ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1703ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1704ccaff030SJeremy L Thompson 
1705ccaff030SJeremy L Thompson   // Clean up PETSc
1706ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1707ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1708ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1709ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1710ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1711ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1712ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1713ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1714ccaff030SJeremy L Thompson   return PetscFinalize();
1715ccaff030SJeremy L Thompson }
1716ccaff030SJeremy L Thompson 
1717