xref: /libCEED/examples/fluids/navierstokes.c (revision 9d801c56255f264a51a1fab570ac18297502b3d3)
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,
1568b982baeSLeila Ghaffari     .applySur_rhs              = DC_Sur,
1578b982baeSLeila 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;
222cfa64770SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs_sur[6], op_rhs,
223cfa64770SLeila Ghaffari                op_ifunction_vol, op_ifunction_sur[6], op_ifunction;
224ccaff030SJeremy L Thompson   Vec M;
225ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
226ccaff030SJeremy L Thompson   PetscInt contsteps;
227ccaff030SJeremy L Thompson };
228ccaff030SJeremy L Thompson 
229ccaff030SJeremy L Thompson struct Units_ {
230ccaff030SJeremy L Thompson   // fundamental units
231ccaff030SJeremy L Thompson   PetscScalar meter;
232ccaff030SJeremy L Thompson   PetscScalar kilogram;
233ccaff030SJeremy L Thompson   PetscScalar second;
234ccaff030SJeremy L Thompson   PetscScalar Kelvin;
235ccaff030SJeremy L Thompson   // derived units
236ccaff030SJeremy L Thompson   PetscScalar Pascal;
237ccaff030SJeremy L Thompson   PetscScalar JperkgK;
238ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
239ccaff030SJeremy L Thompson   PetscScalar WpermK;
240ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
241ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
242ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
243ccaff030SJeremy L Thompson };
244ccaff030SJeremy L Thompson 
245ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
246ccaff030SJeremy L Thompson struct SimpleBC_ {
247cfa64770SLeila Ghaffari   PetscInt nwall, nslip[3], noutflow;
24884d34d69SLeila Ghaffari   PetscInt walls[6], slips[3][6], outflow[6];
24984d34d69SLeila Ghaffari   PetscBool userbc;
250ccaff030SJeremy L Thompson };
251ccaff030SJeremy L Thompson 
252ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
253ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
254ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
255ccaff030SJeremy L Thompson }
256ccaff030SJeremy L Thompson 
257ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
258ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
25984d34d69SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
260ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
261ccaff030SJeremy L Thompson 
262ccaff030SJeremy L Thompson   PetscSection   section;
2630c6c0b13SLeila Ghaffari   PetscInt       p, Nelem, Ndof, *erestrict, eoffset, nfields, dim,
2640c6c0b13SLeila Ghaffari                  depth;
2650c6c0b13SLeila Ghaffari   DMLabel        depthLabel;
2660c6c0b13SLeila Ghaffari   IS             depthIS, iterIS;
26784d34d69SLeila Ghaffari   Vec            Uloc;
2680c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
269ccaff030SJeremy L Thompson   PetscErrorCode ierr;
270ccaff030SJeremy L Thompson 
271ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
272ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
273da51bdd9SLeila Ghaffari   dim -= height;
274ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
275ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
276ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
277ccaff030SJeremy L Thompson   fieldoff[0] = 0;
278ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
279ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
280ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
281ccaff030SJeremy L Thompson   }
282ccaff030SJeremy L Thompson 
2830c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2840c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2850c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2860c6c0b13SLeila Ghaffari   if (domainLabel) {
2870c6c0b13SLeila Ghaffari     IS domainIS;
2880c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2890c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2900c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2910c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2920c6c0b13SLeila Ghaffari   } else {
2930c6c0b13SLeila Ghaffari     iterIS = depthIS;
2940c6c0b13SLeila Ghaffari   }
2950c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2960c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
297ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2980c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2990c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
300ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
30184d34d69SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
30284d34d69SLeila Ghaffari                                    &numindices, &indices, NULL, NULL);
30384d34d69SLeila Ghaffari     CHKERRQ(ierr);
30484d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
30584d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
30684d34d69SLeila Ghaffari           c);
307ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
308ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
309ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
310ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
311ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
312ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
313ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
314ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
315ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
316ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
317ccaff030SJeremy L Thompson                      c, i, f, j);
318ccaff030SJeremy L Thompson         }
319ccaff030SJeremy L Thompson       }
320ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
321ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
3226f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
323ccaff030SJeremy L Thompson     }
32484d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
32584d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
32684d34d69SLeila Ghaffari     CHKERRQ(ierr);
327ccaff030SJeremy L Thompson   }
3280c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3290c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3300c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
331ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3320c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3330c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3340c6c0b13SLeila Ghaffari 
335ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
336ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
337ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3386f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3396f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3406f55dfd5Svaleriabarra                             Erestrict);
341ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
342ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
343ccaff030SJeremy L Thompson }
344ccaff030SJeremy L Thompson 
345c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
346bd910870SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt ncompx, CeedInt dim,
347cfa64770SLeila Ghaffari     CeedInt height, DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q,
348c96c872fSLeila Ghaffari     CeedInt qdatasize, CeedElemRestriction *restrictq,
349c96c872fSLeila Ghaffari     CeedElemRestriction *restrictx, CeedElemRestriction *restrictqdi) {
350c96c872fSLeila Ghaffari 
351c96c872fSLeila Ghaffari   DM dmcoord;
352c96c872fSLeila Ghaffari   CeedInt localNelem;
353c96c872fSLeila Ghaffari   CeedInt Qdim = CeedIntPow(Q, dim);
354c96c872fSLeila Ghaffari   PetscErrorCode ierr;
355c96c872fSLeila Ghaffari 
356c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
357c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
358c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
359c96c872fSLeila Ghaffari   CHKERRQ(ierr);
360c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
361c96c872fSLeila Ghaffari   CHKERRQ(ierr);
362c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
363c96c872fSLeila Ghaffari   CHKERRQ(ierr);
364c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
365c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
366c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
367c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
368c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
369c96c872fSLeila Ghaffari }
370c96c872fSLeila Ghaffari 
371ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
372ccaff030SJeremy L Thompson   PetscErrorCode ierr;
373ccaff030SJeremy L Thompson   PetscInt m;
374ccaff030SJeremy L Thompson 
375ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
376ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
377ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
378ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
379ccaff030SJeremy L Thompson }
380ccaff030SJeremy L Thompson 
381ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
382ccaff030SJeremy L Thompson   PetscErrorCode ierr;
383ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
384ccaff030SJeremy L Thompson   PetscScalar *a;
385ccaff030SJeremy L Thompson 
386ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
387ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
388ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
389ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
39084d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
39184d34d69SLeila Ghaffari                                   mpetsc, mceed);
392ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
393ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
394ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
395ccaff030SJeremy L Thompson }
396ccaff030SJeremy L Thompson 
397ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
398ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
399ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
400ccaff030SJeremy L Thompson   PetscErrorCode ierr;
401ccaff030SJeremy L Thompson   Vec Qbc;
402ccaff030SJeremy L Thompson 
403ccaff030SJeremy L Thompson   PetscFunctionBegin;
404ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
405ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
406ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
408ccaff030SJeremy L Thompson }
409ccaff030SJeremy L Thompson 
410ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
411ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
412ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
413ccaff030SJeremy L Thompson   PetscErrorCode ierr;
414ccaff030SJeremy L Thompson   User user = *(User *)userData;
415ccaff030SJeremy L Thompson   PetscScalar *q, *g;
416ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
417ccaff030SJeremy L Thompson 
418ccaff030SJeremy L Thompson   // Global-to-local
419ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
420ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
421ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
422ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
423ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
424ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
425ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
426ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
427ccaff030SJeremy L Thompson 
428ccaff030SJeremy L Thompson   // Ceed Vectors
429ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
430ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
431ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
432ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
433ccaff030SJeremy L Thompson 
434ccaff030SJeremy L Thompson   // Apply CEED operator
435ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
436ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
437ccaff030SJeremy L Thompson 
438ccaff030SJeremy L Thompson   // Restore vectors
439ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
440ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson 
442ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
443ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
444ccaff030SJeremy L Thompson 
445ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
446ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
447ccaff030SJeremy L Thompson   CHKERRQ(ierr);
448ccaff030SJeremy L Thompson 
449ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
450ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
451ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
452ccaff030SJeremy L Thompson }
453ccaff030SJeremy L Thompson 
454ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
455ccaff030SJeremy L Thompson                                    void *userData) {
456ccaff030SJeremy L Thompson   PetscErrorCode ierr;
457ccaff030SJeremy L Thompson   User user = *(User *)userData;
458ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
459ccaff030SJeremy L Thompson   PetscScalar *g;
460ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
461ccaff030SJeremy L Thompson 
462ccaff030SJeremy L Thompson   // Global-to-local
463ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
464ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
465ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
466ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
468ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
469ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
470ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
471ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
472ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
473ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
474ccaff030SJeremy L Thompson 
475ccaff030SJeremy L Thompson   // Ceed Vectors
476ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
479ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
480ccaff030SJeremy L Thompson                      (PetscScalar *)q);
481ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
482ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
483ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
484ccaff030SJeremy L Thompson 
485ccaff030SJeremy L Thompson   // Apply CEED operator
486ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
487ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
488ccaff030SJeremy L Thompson 
489ccaff030SJeremy L Thompson   // Restore vectors
490ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
491ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson 
494ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
495ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
496ccaff030SJeremy L Thompson 
497ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
499ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
500ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
501ccaff030SJeremy L Thompson }
502ccaff030SJeremy L Thompson 
503ccaff030SJeremy L Thompson // User provided TS Monitor
504ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
505ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
506ccaff030SJeremy L Thompson   User user = ctx;
507ccaff030SJeremy L Thompson   Vec Qloc;
508ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
509ccaff030SJeremy L Thompson   PetscViewer viewer;
510ccaff030SJeremy L Thompson   PetscErrorCode ierr;
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson   // Set up output
513ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
514ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
515ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
516ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
517ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
518ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
520ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
521ccaff030SJeremy L Thompson 
522ccaff030SJeremy L Thompson   // Output
523ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
524ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
525ccaff030SJeremy L Thompson   CHKERRQ(ierr);
526ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
527ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
528ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
529*9d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
530ccaff030SJeremy L Thompson   if (user->dmviz) {
531ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
532ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
533ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
534ccaff030SJeremy L Thompson 
535ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
536ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
538ccaff030SJeremy L Thompson     CHKERRQ(ierr);
539ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
542ccaff030SJeremy L Thompson     CHKERRQ(ierr);
543ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
544ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
545ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
546ccaff030SJeremy L Thompson     CHKERRQ(ierr);
547ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
548ccaff030SJeremy L Thompson                               filepath_refined,
549ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
550ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
553ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
554ccaff030SJeremy L Thompson   }
555ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
556ccaff030SJeremy L Thompson 
557ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
558ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
559ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
560ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
561ccaff030SJeremy L Thompson   CHKERRQ(ierr);
562ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
563ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson 
565ccaff030SJeremy L Thompson   // Save time stamp
566ccaff030SJeremy L Thompson   // Dimensionalize time back
567ccaff030SJeremy L Thompson   time /= user->units->second;
568ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
569ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
570ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
571ccaff030SJeremy L Thompson   CHKERRQ(ierr);
572ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
573ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
574ccaff030SJeremy L Thompson   #else
575ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
576ccaff030SJeremy L Thompson   #endif
577ccaff030SJeremy L Thompson   CHKERRQ(ierr);
578ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
579ccaff030SJeremy L Thompson 
580ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
581ccaff030SJeremy L Thompson }
582ccaff030SJeremy L Thompson 
58384d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
584ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
585ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
586ccaff030SJeremy L Thompson   PetscErrorCode ierr;
587ccaff030SJeremy L Thompson   CeedVector multlvec;
588ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
589ccaff030SJeremy L Thompson 
590ccaff030SJeremy L Thompson   ctxSetup->time = time;
591ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
593ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
594ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
595ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
596ccaff030SJeremy L Thompson 
597ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
598ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
599ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
600ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
601ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
602ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
603ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
604ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
605ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
606ccaff030SJeremy L Thompson   CHKERRQ(ierr);
607ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
609ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
610ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
613ccaff030SJeremy L Thompson }
614ccaff030SJeremy L Thompson 
615ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
616ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
617ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
618ccaff030SJeremy L Thompson   PetscErrorCode ierr;
619ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
620ccaff030SJeremy L Thompson   CeedOperator op_mass;
621ccaff030SJeremy L Thompson   CeedVector mceed;
622ccaff030SJeremy L Thompson   Vec Mloc;
623ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
624ccaff030SJeremy L Thompson 
625ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
626ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
627ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
628ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
629ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
630ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
631ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
632ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
633ccaff030SJeremy L Thompson 
634ccaff030SJeremy L Thompson   // Create the mass operator
635ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
636ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
637ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
638ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
639ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
640ccaff030SJeremy L Thompson 
641ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
643ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
644ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson 
646ccaff030SJeremy L Thompson   {
647ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
648ccaff030SJeremy L Thompson     CeedVector onesvec;
649ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
650ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
651ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
652ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
653ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
654ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
655ccaff030SJeremy L Thompson   }
656ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
657ccaff030SJeremy L Thompson 
658ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
659ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
660ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
661ccaff030SJeremy L Thompson 
662ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
663ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
664ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
665ccaff030SJeremy L Thompson }
666ccaff030SJeremy L Thompson 
66784d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
668ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
669ccaff030SJeremy L Thompson   PetscErrorCode ierr;
670ccaff030SJeremy L Thompson 
671ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
672ccaff030SJeremy L Thompson   {
673ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
674ccaff030SJeremy L Thompson     PetscFE fe;
675ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
676ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
677ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
67832ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
679ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
680ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
681ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
68207af6069Svaleriabarra     {
68307af6069Svaleriabarra       PetscInt comps[1] = {1};
68407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
68507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
68607af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
68707af6069Svaleriabarra       comps[0] = 2;
68807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
68907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
69007af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
69107af6069Svaleriabarra       comps[0] = 3;
69207af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
69307af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
69407af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
69507af6069Svaleriabarra     }
69684d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
69784d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
69884d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
69984d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
70084d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
70184d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
70284d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
70384d34d69SLeila Ghaffari 
70484d34d69SLeila Ghaffari           }
70584d34d69SLeila Ghaffari         }
70684d34d69SLeila Ghaffari       }
70784d34d69SLeila Ghaffari     }
70884d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
70984d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
71084d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
71184d34d69SLeila Ghaffari     {
71284d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
71384d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
71484d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
71584d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
71684d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
71784d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
71884d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
71984d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
72084d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
72184d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
72284d34d69SLeila Ghaffari       } else
72384d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
72484d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
72584d34d69SLeila Ghaffari     }
726ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
727ccaff030SJeremy L Thompson     CHKERRQ(ierr);
728ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
729ccaff030SJeremy L Thompson   }
730ccaff030SJeremy L Thompson   {
731ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
732ccaff030SJeremy L Thompson     PetscSection section;
733ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
734ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
735ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
736ccaff030SJeremy L Thompson     CHKERRQ(ierr);
737ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
738ccaff030SJeremy L Thompson     CHKERRQ(ierr);
739ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
740ccaff030SJeremy L Thompson     CHKERRQ(ierr);
741ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
742ccaff030SJeremy L Thompson     CHKERRQ(ierr);
743ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
744ccaff030SJeremy L Thompson     CHKERRQ(ierr);
745ccaff030SJeremy L Thompson   }
746ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
747ccaff030SJeremy L Thompson }
748ccaff030SJeremy L Thompson 
749ccaff030SJeremy L Thompson int main(int argc, char **argv) {
750ccaff030SJeremy L Thompson   PetscInt ierr;
751ccaff030SJeremy L Thompson   MPI_Comm comm;
75284d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
753ccaff030SJeremy L Thompson   Mat interpviz;
754ccaff030SJeremy L Thompson   TS ts;
755ccaff030SJeremy L Thompson   TSAdapt adapt;
756ccaff030SJeremy L Thompson   User user;
757ccaff030SJeremy L Thompson   Units units;
758ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
75984d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
760ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
761ccaff030SJeremy L Thompson   PetscMPIInt rank;
762ccaff030SJeremy L Thompson   PetscScalar ftime;
763ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
764ccaff030SJeremy L Thompson   Ceed ceed;
76584d34d69SLeila Ghaffari   CeedInt numP, numQ;
766cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
76784d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
76884d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
769cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
770cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
771ccaff030SJeremy L Thompson   CeedScalar Rd;
77284d34d69SLeila Ghaffari   CeedMemType memtyperequested;
773ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
774ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
775ccaff030SJeremy L Thompson   problemType problemChoice;
776ccaff030SJeremy L Thompson   problemData *problem = NULL;
777ccaff030SJeremy L Thompson   StabilizationType stab;
77884d34d69SLeila Ghaffari   testType testChoice;
77984d34d69SLeila Ghaffari   testData *test = NULL;
78084d34d69SLeila Ghaffari   PetscBool implicit;
781cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
782ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
78384d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
78484d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
785ccaff030SJeremy L Thompson   };
786ccaff030SJeremy L Thompson   double start, cpu_time_used;
78784d34d69SLeila Ghaffari   // Check PETSc CUDA support
78884d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
78984d34d69SLeila Ghaffari   // *INDENT-OFF*
79084d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
79184d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
79284d34d69SLeila Ghaffari   #else
79384d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
79484d34d69SLeila Ghaffari   #endif
79584d34d69SLeila Ghaffari   // *INDENT-ON*
796ccaff030SJeremy L Thompson 
797ccaff030SJeremy L Thompson   // Create the libCEED contexts
798ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
799ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
800ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
801ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
802ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
803ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
804ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
805ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
806ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
807ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
808ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
809ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
810ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
811ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
812ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
813ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
814ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
815ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
816ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
817ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
818ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
819ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
820ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
821ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
822ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
823ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
82484d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
82584d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
826ea6e0f84SLeila Ghaffari   PetscInt degreeSur     = 1;        // -
827ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
828ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
829ccaff030SJeremy L Thompson 
830ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
831ccaff030SJeremy L Thompson   if (ierr) return ierr;
832ccaff030SJeremy L Thompson 
833ccaff030SJeremy L Thompson   // Allocate PETSc context
834ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
835ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
836ccaff030SJeremy L Thompson 
837ccaff030SJeremy L Thompson   // Parse command line options
838ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
839ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
840ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
841ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
842ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
843ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
84484d34d69SLeila Ghaffari   testChoice = TEST_NONE;
84584d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
84684d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
84784d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
84884d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
84984d34d69SLeila Ghaffari   test = &testOptions[testChoice];
850ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
851ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
852ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
853ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
854ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
855ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
856ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
857ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
858ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
859ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
860ccaff030SJeremy L Thompson   CHKERRQ(ierr);
86184d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
86284d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
86384d34d69SLeila Ghaffari     CHKERRQ(ierr);
86484d34d69SLeila Ghaffari   }
865ccaff030SJeremy L Thompson   {
866cfa64770SLeila Ghaffari     PetscInt len, len1;
867cfa64770SLeila Ghaffari     PetscBool flg, flg1;
868ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
869ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
870ccaff030SJeremy L Thompson                                 NULL, bc.walls,
871ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
872ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
873ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
874ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
875ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
876ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
877ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
878ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
879ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
880ccaff030SJeremy L Thompson                                    &len), &flg);
881ccaff030SJeremy L Thompson       CHKERRQ(ierr);
88284d34d69SLeila Ghaffari       if (flg) {
88384d34d69SLeila Ghaffari         bc.nslip[j] = len;
88484d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
88584d34d69SLeila Ghaffari       }
886ccaff030SJeremy L Thompson     }
887cfa64770SLeila Ghaffari     ierr = PetscOptionsIntArray("-bc_outflow",
888cfa64770SLeila Ghaffari                               "Use outflow boundary conditions on this list of faces",
889cfa64770SLeila Ghaffari                               NULL, bc.outflow,
890cfa64770SLeila Ghaffari                               (len1 = sizeof(bc.outflow) / sizeof(bc.outflow[0]),
891cfa64770SLeila Ghaffari                               &len1), &flg1); CHKERRQ(ierr);
892cfa64770SLeila Ghaffari     if (flg1) bc.noutflow = len1;
893ccaff030SJeremy L Thompson   }
894cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
895cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
896cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
897ccaff030SJeremy L Thompson   CHKERRQ(ierr);
898ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
899ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
900ccaff030SJeremy L Thompson   meter = fabs(meter);
901ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
902ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
903ccaff030SJeremy L Thompson   second = fabs(second);
904ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
905ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
906ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
907ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
908ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
909ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
910ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
911ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
912ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
913ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
914ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
915ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
916ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
917ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
918ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
919ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
920ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
921ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
922ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
923ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
924ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
925ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
926ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
927ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
928ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
929ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
930ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
931ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
932ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
933ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
934ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
93584d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
93684d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
93784d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
93884d34d69SLeila Ghaffari     CHKERRQ(ierr);
93984d34d69SLeila Ghaffari   }
940ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
941ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
942ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
943ccaff030SJeremy L Thompson   CHKERRQ(ierr);
94484d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
94584d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
94684d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
94784d34d69SLeila Ghaffari     CHKERRQ(ierr);
94884d34d69SLeila Ghaffari   }
949ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
950ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
951ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
952ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
953ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
954ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
955ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
956ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
957ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
958ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
959ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
960ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
961ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
962ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
963ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
964ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
965ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
966ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
967ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
968ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
969ccaff030SJeremy L Thompson   n = problem->dim;
970ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
971ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
972ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
973ccaff030SJeremy L Thompson   {
974ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
975ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
976ccaff030SJeremy L Thompson     if (norm > 0) {
977ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
978ccaff030SJeremy L Thompson     }
979ccaff030SJeremy L Thompson   }
980ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
981ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
982ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
983ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
984ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
98584d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
98684d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
98784d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
98884d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
98984d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
990ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
991ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
992ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
99384d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
99484d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
99584d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
99684d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
99784d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
99884d34d69SLeila Ghaffari   CHKERRQ(ierr);
999ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1000ccaff030SJeremy L Thompson 
1001ccaff030SJeremy L Thompson   // Define derived units
1002ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1003ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1004ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1005ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1006ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1007ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1008ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
1009ccaff030SJeremy L Thompson 
1010ccaff030SJeremy L Thompson   // Scale variables to desired units
1011ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1012ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1013ccaff030SJeremy L Thompson   P0 *= Pascal;
1014ccaff030SJeremy L Thompson   N *= (1./second);
1015ccaff030SJeremy L Thompson   cv *= JperkgK;
1016ccaff030SJeremy L Thompson   cp *= JperkgK;
1017ccaff030SJeremy L Thompson   Rd = cp - cv;
1018ccaff030SJeremy L Thompson   g *= mpersquareds;
1019ccaff030SJeremy L Thompson   mu *= Pascal * second;
1020ccaff030SJeremy L Thompson   k *= WpermK;
1021ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1022ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1023ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1024ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1025ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1026ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1027ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1028ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1029ccaff030SJeremy L Thompson 
1030ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1031cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1032ccaff030SJeremy L Thompson   // Set up the libCEED context
1033ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1034ccaff030SJeremy L Thompson     .theta0 = theta0,
1035ccaff030SJeremy L Thompson     .thetaC = thetaC,
1036ccaff030SJeremy L Thompson     .P0 = P0,
1037ccaff030SJeremy L Thompson     .N = N,
1038ccaff030SJeremy L Thompson     .cv = cv,
1039ccaff030SJeremy L Thompson     .cp = cp,
1040ccaff030SJeremy L Thompson     .Rd = Rd,
1041ccaff030SJeremy L Thompson     .g = g,
1042ccaff030SJeremy L Thompson     .rc = rc,
1043ccaff030SJeremy L Thompson     .lx = lx,
1044ccaff030SJeremy L Thompson     .ly = ly,
1045ccaff030SJeremy L Thompson     .lz = lz,
1046ccaff030SJeremy L Thompson     .center[0] = center[0],
1047ccaff030SJeremy L Thompson     .center[1] = center[1],
1048ccaff030SJeremy L Thompson     .center[2] = center[2],
1049ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1050ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1051ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
1052ccaff030SJeremy L Thompson     .time = 0,
1053ccaff030SJeremy L Thompson   };
1054ccaff030SJeremy L Thompson 
105584d34d69SLeila Ghaffari   // Create the mesh
1056ccaff030SJeremy L Thompson   {
1057ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1058ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
105984d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1060ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1061ccaff030SJeremy L Thompson   }
106284d34d69SLeila Ghaffari 
106384d34d69SLeila Ghaffari   // Distribute the mesh over processes
106484d34d69SLeila Ghaffari   {
1065ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1066ccaff030SJeremy L Thompson     PetscPartitioner part;
1067ccaff030SJeremy L Thompson 
1068ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1069ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1070ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1071ccaff030SJeremy L Thompson     if (dmDist) {
1072ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1073ccaff030SJeremy L Thompson       dm  = dmDist;
1074ccaff030SJeremy L Thompson     }
1075ccaff030SJeremy L Thompson   }
1076ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1077ccaff030SJeremy L Thompson 
107884d34d69SLeila Ghaffari   // Setup DM
1079ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1080ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
108184d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
108284d34d69SLeila Ghaffari 
108384d34d69SLeila Ghaffari   // Refine DM for high-order viz
1084ccaff030SJeremy L Thompson   dmviz = NULL;
1085ccaff030SJeremy L Thompson   interpviz = NULL;
1086ccaff030SJeremy L Thompson   if (viz_refine) {
1087ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1088ff6701fcSJed Brown 
1089ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1090ff6701fcSJed Brown     dmhierarchy[0] = dm;
109184d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1092ff6701fcSJed Brown       Mat interp_next;
1093ff6701fcSJed Brown 
1094ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1095ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1096ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1097ff6701fcSJed Brown       d = (d + 1) / 2;
1098ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1099ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1100ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1101ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1102ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1103ff6701fcSJed Brown       else {
1104ff6701fcSJed Brown         Mat C;
1105ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1106ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1107ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1108ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1109ff6701fcSJed Brown         interpviz = C;
1110ff6701fcSJed Brown       }
1111ff6701fcSJed Brown     }
1112cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1113ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1114cb3e2689Svaleriabarra     }
1115ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1116ccaff030SJeremy L Thompson   }
1117ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1118ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1119ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1120ccaff030SJeremy L Thompson   lnodes /= ncompq;
1121ccaff030SJeremy L Thompson 
112284d34d69SLeila Ghaffari   // Initialize CEED
112384d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
112484d34d69SLeila Ghaffari   // Set memtype
112584d34d69SLeila Ghaffari   CeedMemType memtypebackend;
112684d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
112784d34d69SLeila Ghaffari   // Check memtype compatibility
112884d34d69SLeila Ghaffari   if (!setmemtyperequest)
112984d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
113084d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
113184d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
113284d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
113384d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
113484d34d69SLeila Ghaffari 
113584d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
113684d34d69SLeila Ghaffari   numP = degree + 1;
113784d34d69SLeila Ghaffari   numQ = numP + qextra;
113884d34d69SLeila Ghaffari 
113984d34d69SLeila Ghaffari     // Print summary
114084d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1141ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1142ccaff030SJeremy L Thompson     int comm_size;
1143ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1144ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1145ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
114684d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1147ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1148ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1149ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
115084d34d69SLeila Ghaffari     const char *usedresource;
115184d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1152ccaff030SJeremy L Thompson 
115384d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
115484d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
115584d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
115684d34d69SLeila Ghaffari                        "  Problem:\n"
115784d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
115884d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
115984d34d69SLeila Ghaffari                        "  PETSc:\n"
116084d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
116184d34d69SLeila Ghaffari                        "  libCEED:\n"
116284d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
116384d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
116484d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
116584d34d69SLeila Ghaffari                        "  Mesh:\n"
116684d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
116784d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
116884d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
116984d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
117084d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
117184d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
117284d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
117384d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
117484d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
117584d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
117684d34d69SLeila Ghaffari                        (setmemtyperequest) ?
117784d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
117884d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
117984d34d69SLeila Ghaffari     CHKERRQ(ierr);
11800c6c0b13SLeila Ghaffari   }
11810c6c0b13SLeila Ghaffari 
1182ccaff030SJeremy L Thompson   // Set up global mass vector
1183ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1184ccaff030SJeremy L Thompson 
118584d34d69SLeila Ghaffari   // Set up libCEED
1186ccaff030SJeremy L Thompson   // CEED Bases
1187ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
118884d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
118984d34d69SLeila Ghaffari                                   &basisq);
119084d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
119184d34d69SLeila Ghaffari                                   &basisx);
119284d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
119384d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1194ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1195ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1196ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1197ccaff030SJeremy L Thompson 
1198ccaff030SJeremy L Thompson   // CEED Restrictions
119984d34d69SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dim, 0, 0, 0, numP, numQ,
120084d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
120184d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1202ccaff030SJeremy L Thompson 
1203ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1204ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1205ccaff030SJeremy L Thompson 
1206ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1207bd910870SLeila Ghaffari   CeedInt NqptsVol;
120884d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
120984d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
12108b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
121184d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1212ccaff030SJeremy L Thompson 
1213ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1214ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1215ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1216ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1217ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
12188b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1219ccaff030SJeremy L Thompson 
1220ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1221ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1222ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1223ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1224ccaff030SJeremy L Thompson 
1225ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1226ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1227ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1228ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1229ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1230ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
12318b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1232ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1233ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1234ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1235ccaff030SJeremy L Thompson   }
1236ccaff030SJeremy L Thompson 
1237ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1238ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1239ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1240ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1241ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1242ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1243ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
12448b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1245ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1246ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1247ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1248ccaff030SJeremy L Thompson   }
1249ccaff030SJeremy L Thompson 
1250ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1251ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
125284d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1253ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
125484d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
125584d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1256ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1257ccaff030SJeremy L Thompson 
1258ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1259ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
126084d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
126184d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1262ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1263ccaff030SJeremy L Thompson 
126484d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
126584d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
126684d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1267ccaff030SJeremy L Thompson 
1268ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1269ccaff030SJeremy L Thompson     CeedOperator op;
1270ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
127184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
127284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
127384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
12748b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
127584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
127684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
127784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1278ccaff030SJeremy L Thompson     user->op_rhs = op;
1279ccaff030SJeremy L Thompson   }
1280ccaff030SJeremy L Thompson 
1281ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1282ccaff030SJeremy L Thompson     CeedOperator op;
1283ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
128484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
128584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
128684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
128784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
12888b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
128984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
129084d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
129184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1292ccaff030SJeremy L Thompson     user->op_ifunction = op;
1293ccaff030SJeremy L Thompson   }
1294ccaff030SJeremy L Thompson 
1295cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1296cfa64770SLeila Ghaffari   // Outflow Boundary Condition
12976a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
12986a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
1299cfa64770SLeila Ghaffari   CeedInt numP_Sur, numQ_Sur;
13006a0edaf9SLeila Ghaffari   CeedInt height = 1;
13016a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
13026a0edaf9SLeila Ghaffari   numP_Sur = degreeSur + 1;
13036a0edaf9SLeila Ghaffari   numQ_Sur = numP_Sur + qextraSur;
1304cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1305cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1306cfa64770SLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
1307cfa64770SLeila Ghaffari   CeedInt NqptsSur;
1308cfa64770SLeila Ghaffari   PetscInt localNelemSur[6];
1309cfa64770SLeila Ghaffari   CeedVector qdataSur[6], qdataSur_[6];
1310cfa64770SLeila Ghaffari   CeedQFunction qf_setupSur, qf_rhsSur, qf_ifunctionSur;
1311cfa64770SLeila Ghaffari   CeedOperator op_setupSur[6], op_setupSur_[6];
1312cfa64770SLeila Ghaffari   PetscInt numOutFlow = bc.noutflow;
1313cfa64770SLeila Ghaffari   DMLabel domainLabel;
1314cfa64770SLeila Ghaffari 
1315cfa64770SLeila Ghaffari   // Get Label for the boundaries
1316cfa64770SLeila Ghaffari   ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
1317cfa64770SLeila Ghaffari 
1318cfa64770SLeila Ghaffari   // CEED bases for the boundaries
13196a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
13206a0edaf9SLeila Ghaffari                                   &basisqSur);
13216a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
13226a0edaf9SLeila Ghaffari                                   &basisxSur);
13236a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
13246a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
13256a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
13266a0edaf9SLeila Ghaffari 
1327cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
13286a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
13296a0edaf9SLeila Ghaffari                               &qf_setupSur);
13306a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
13316a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
13326a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
13336a0edaf9SLeila Ghaffari 
13346a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
1335cfa64770SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator on the Surface
13366a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
13376a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
13386a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
13398b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "dq", ncompq*dim, CEED_EVAL_GRAD); //assumed volumetric elements
13406a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
13416a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
1342cfa64770SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdata", qdatasizeVol, CEED_EVAL_NONE); //assumed volumetric elements
13436a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
13446a0edaf9SLeila Ghaffari   }
13456a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
13466a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
13476a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
13486a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
13496a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
13508b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "dq", ncompq*dim, CEED_EVAL_GRAD); //assumed volumetric elements
13516a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
13526a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
1353cfa64770SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdata", qdatasizeVol, CEED_EVAL_NONE); //assumed volumetric elements
13546a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
13556a0edaf9SLeila Ghaffari   }
1356cfa64770SLeila Ghaffari   // Create CEED Operator for each face
1357cfa64770SLeila Ghaffari   for(CeedInt i=0; i<numOutFlow; i++){
1358cfa64770SLeila Ghaffari     ierr = GetRestrictionForDomain(ceed, dm, ncompx, dimSur, height, domainLabel, bc.outflow[i], numP_Sur,
1359cfa64770SLeila Ghaffari                                    numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
1360cfa64770SLeila Ghaffari                                    &restrictqdiSur[i]); CHKERRQ(ierr);
1361cfa64770SLeila Ghaffari     // Create the CEED vectors that will be needed in setup
1362cfa64770SLeila Ghaffari     CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
1363cfa64770SLeila Ghaffari     CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
1364cfa64770SLeila Ghaffari     CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur_[i]);
13656a0edaf9SLeila Ghaffari 
13666a0edaf9SLeila Ghaffari     // Create the operator that builds the quadrature data for the NS operator
1367cfa64770SLeila Ghaffari     CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
1368cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, CEED_VECTOR_ACTIVE);
1369cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
13706a0edaf9SLeila Ghaffari                          basisxSur, CEED_VECTOR_NONE);
1371cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
13726a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
13736a0edaf9SLeila Ghaffari 
1374cfa64770SLeila Ghaffari     // Utility operator that builds the quadrature data for computing viscous terms
1375cfa64770SLeila Ghaffari     CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupSur_[i]);
137684d34d69SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "dx", restrictxSur[i], basisx, CEED_VECTOR_ACTIVE);
1377cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "weight", CEED_ELEMRESTRICTION_NONE,
137884d34d69SLeila Ghaffari                          basisx, CEED_VECTOR_NONE);
1379cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "qdata", restrictqdiSur[i],
1380cfa64770SLeila Ghaffari                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
13816a0edaf9SLeila Ghaffari 
13826a0edaf9SLeila Ghaffari     if (qf_rhsSur) { // Create the RHS physics operator
13836a0edaf9SLeila Ghaffari       CeedOperator op;
13846a0edaf9SLeila Ghaffari       CeedOperatorCreate(ceed, qf_rhsSur, NULL, NULL, &op);
1385cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1386cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "dq", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1387cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdataSur", restrictqdiSur[i],
1388cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
1389cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "x", restrictxSur[i], basisxSur, xcorners);
1390cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdata", restrictqdiSur[i],
1391cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur_[i]);
1392cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1393cfa64770SLeila Ghaffari       user->op_rhs_sur[i] = op;
13946a0edaf9SLeila Ghaffari     }
13956a0edaf9SLeila Ghaffari 
13966a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) { // Create the IFunction operator
13976a0edaf9SLeila Ghaffari       CeedOperator op;
13986a0edaf9SLeila Ghaffari       CeedOperatorCreate(ceed, qf_ifunctionSur, NULL, NULL, &op);
1399cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1400cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "dq", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1401cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdataSur", restrictqdiSur[i],
1402cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
1403cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "x", restrictxSur[i], basisxSur, xcorners);
1404cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdata", restrictqdiSur[i],
1405cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur_[i]);
1406cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1407cfa64770SLeila Ghaffari       user->op_ifunction_sur[i] = op;
1408cfa64770SLeila Ghaffari     }
14096a0edaf9SLeila Ghaffari   }
14106a0edaf9SLeila Ghaffari   // Composite Operaters
141184d34d69SLeila Ghaffari   // IFunction
14126a0edaf9SLeila Ghaffari   if (user->op_ifunction_vol) {
1413cfa64770SLeila Ghaffari     if (numOutFlow>0) {
14146a0edaf9SLeila Ghaffari       // Composite Operators for the IFunction
14156a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_ifunction);
14166a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol);
1417cfa64770SLeila Ghaffari       for(CeedInt i=0; i<numOutFlow; i++){
1418cfa64770SLeila Ghaffari         CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_sur[i]);
1419cfa64770SLeila Ghaffari       }
14206a0edaf9SLeila Ghaffari   } else {
14216a0edaf9SLeila Ghaffari     user->op_ifunction = user->op_ifunction_vol;
14226a0edaf9SLeila Ghaffari     }
14236a0edaf9SLeila Ghaffari   }
142484d34d69SLeila Ghaffari   // RHS
14256a0edaf9SLeila Ghaffari   if (user->op_rhs_vol) {
1426cfa64770SLeila Ghaffari     if (numOutFlow == 1) {
14276a0edaf9SLeila Ghaffari       // Composite Operators for the RHS
14286a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_rhs);
14296a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_vol);
1430cfa64770SLeila Ghaffari       for(CeedInt i=0; i<numOutFlow; i++){
1431cfa64770SLeila Ghaffari         CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_sur[i]);
1432cfa64770SLeila Ghaffari       }
14336a0edaf9SLeila Ghaffari   } else {
14346a0edaf9SLeila Ghaffari     user->op_rhs = user->op_rhs_vol;
14356a0edaf9SLeila Ghaffari     }
14366a0edaf9SLeila Ghaffari   }
1437cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1438ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1439ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1440ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1441ccaff030SJeremy L Thompson     .CtauS = CtauS,
1442ccaff030SJeremy L Thompson     .strong_form = strong_form,
1443ccaff030SJeremy L Thompson     .stabilization = stab,
1444ccaff030SJeremy L Thompson   };
1445ccaff030SJeremy L Thompson   switch (problemChoice) {
1446ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1447ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1448ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1449ccaff030SJeremy L Thompson           sizeof ctxNS);
14506a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
14516a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
14526a0edaf9SLeila Ghaffari           sizeof ctxNS);
1453ccaff030SJeremy L Thompson     break;
1454ccaff030SJeremy L Thompson   case NS_ADVECTION:
1455ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1456ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1457ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1458ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1459ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
14606a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
14616a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
14626a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
14636a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1464ccaff030SJeremy L Thompson   }
1465ccaff030SJeremy L Thompson 
1466ccaff030SJeremy L Thompson   // Set up PETSc context
1467ccaff030SJeremy L Thompson   // Set up units structure
1468ccaff030SJeremy L Thompson   units->meter = meter;
1469ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1470ccaff030SJeremy L Thompson   units->second = second;
1471ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1472ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1473ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1474ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1475ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1476ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1477ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1478ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1479ccaff030SJeremy L Thompson 
1480ccaff030SJeremy L Thompson   // Set up user structure
1481ccaff030SJeremy L Thompson   user->comm = comm;
1482ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1483ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1484ccaff030SJeremy L Thompson   user->units = units;
1485ccaff030SJeremy L Thompson   user->dm = dm;
1486ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1487ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1488ccaff030SJeremy L Thompson   user->ceed = ceed;
1489ccaff030SJeremy L Thompson 
14908b982baeSLeila Ghaffari   // Calculate qdata and ICs
1491ccaff030SJeremy L Thompson   // Set up state global and local vectors
1492ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1493ccaff030SJeremy L Thompson 
1494cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1495ccaff030SJeremy L Thompson 
1496ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1497ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
14988b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
149984d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1500ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1501ccaff030SJeremy L Thompson 
150284d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
150384d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1504ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1505ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1506ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1507ccaff030SJeremy L Thompson     // slower execution.
1508ccaff030SJeremy L Thompson     Vec Qbc;
1509ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1510ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1511ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1512ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1513ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1514ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1515ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
151684d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
151784d34d69SLeila Ghaffari     CHKERRQ(ierr);
1518ccaff030SJeremy L Thompson   }
1519ccaff030SJeremy L Thompson 
1520ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1521ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1522ccaff030SJeremy L Thompson   // Gather initial Q values
1523ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1524ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1525ccaff030SJeremy L Thompson     PetscViewer viewer;
1526ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1527ccaff030SJeremy L Thompson     // Read input
1528ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1529ccaff030SJeremy L Thompson                          user->outputfolder);
1530ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1531ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1532ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1533ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1534ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1535ccaff030SJeremy L Thompson   }
1536ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1537ccaff030SJeremy L Thompson 
1538ccaff030SJeremy L Thompson // Create and setup TS
1539ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1540ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1541ccaff030SJeremy L Thompson   if (implicit) {
1542ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1543ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1544ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1545ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1546ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1547ccaff030SJeremy L Thompson     }
1548ccaff030SJeremy L Thompson   } else {
1549ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1550ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1551ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1552ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1553ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1554ccaff030SJeremy L Thompson   }
1555ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1556ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1557ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
155884d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1559ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1560ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1561ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1562ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1563ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1564ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
156584d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1566ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1567ccaff030SJeremy L Thompson     }
1568ccaff030SJeremy L Thompson   } else { // continue from time of last output
1569ccaff030SJeremy L Thompson     PetscReal time;
1570ccaff030SJeremy L Thompson     PetscInt count;
1571ccaff030SJeremy L Thompson     PetscViewer viewer;
1572ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1573ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1574ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1575ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1576ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1577ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1578ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1579ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1580ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1581ccaff030SJeremy L Thompson   }
158284d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1583ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1584ccaff030SJeremy L Thompson   }
1585ccaff030SJeremy L Thompson 
1586ccaff030SJeremy L Thompson   // Solve
1587ccaff030SJeremy L Thompson   start = MPI_Wtime();
1588ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1589ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1590ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1591ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1592ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1593ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
159484d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1595ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
159684d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1597ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1598ccaff030SJeremy L Thompson   }
1599ccaff030SJeremy L Thompson 
1600ccaff030SJeremy L Thompson   // Get error
160184d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1602ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1603ccaff030SJeremy L Thompson     PetscReal norm;
1604ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1605ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1607ccaff030SJeremy L Thompson 
160884d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
160984d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1610ccaff030SJeremy L Thompson 
1611ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1612ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1613cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1614ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1615ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1616ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
161784d34d69SLeila Ghaffari     // Clean up vectors
161884d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
161984d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1620ccaff030SJeremy L Thompson   }
1621ccaff030SJeremy L Thompson 
1622ccaff030SJeremy L Thompson   // Output Statistics
1623ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
162484d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1625ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1626ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1627ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1628ccaff030SJeremy L Thompson   }
162984d34d69SLeila Ghaffari   // Output numerical values from command line
163084d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
163184d34d69SLeila Ghaffari 
163284d34d69SLeila Ghaffari   // compare reference solution values with current run
163384d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
163484d34d69SLeila Ghaffari     PetscViewer viewer;
163584d34d69SLeila Ghaffari     // Read reference file
163684d34d69SLeila Ghaffari     Vec Qref;
163784d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
163884d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
163984d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
164084d34d69SLeila Ghaffari     CHKERRQ(ierr);
164184d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
164284d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
164384d34d69SLeila Ghaffari 
164484d34d69SLeila Ghaffari     // Compute error with respect to reference solution
164584d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
164684d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
164784d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
164884d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
164984d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
165084d34d69SLeila Ghaffari     // Check error
165184d34d69SLeila Ghaffari     if (error > test->testtol) {
165284d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
165384d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
165484d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
165584d34d69SLeila Ghaffari     }
165684d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
165784d34d69SLeila Ghaffari   }
16589cf88b28Svaleriabarra 
1659ccaff030SJeremy L Thompson   // Clean up libCEED
16608b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1661ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1662ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1663ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1664ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
166584d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
166684d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
166784d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
166884d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
166984d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
167084d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1671ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1672ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1673ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1674ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1675ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1676ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
16776a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
16786a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
16796a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
16806a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
16816a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
16826a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
16836a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
16846a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
16856a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
1686cfa64770SLeila Ghaffari   for(CeedInt i=0; i<numOutFlow; i++){
1687cfa64770SLeila Ghaffari     CeedVectorDestroy(&qdataSur[i]);
1688cfa64770SLeila Ghaffari     CeedVectorDestroy(&qdataSur_[i]);
1689cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictqSur[i]);
1690cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictxSur[i]);
1691cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictqdiSur[i]);
1692cfa64770SLeila Ghaffari     CeedOperatorDestroy(&op_setupSur[i]);
1693cfa64770SLeila Ghaffari     CeedOperatorDestroy(&op_setupSur_[i]);
1694cfa64770SLeila Ghaffari     CeedOperatorDestroy(&user->op_rhs_sur[i]);
1695cfa64770SLeila Ghaffari     CeedOperatorDestroy(&user->op_ifunction_sur[i]);
1696cfa64770SLeila Ghaffari   }
1697ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1698ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1699ccaff030SJeremy L Thompson 
1700ccaff030SJeremy L Thompson   // Clean up PETSc
1701ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1702ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1703ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1704ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1705ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1706ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1707ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1708ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1709ccaff030SJeremy L Thompson   return PetscFinalize();
1710ccaff030SJeremy L Thompson }
1711ccaff030SJeremy L Thompson 
1712