xref: /libCEED/examples/fluids/navierstokes.c (revision 1119eeeee8fcc2d26e311e6320cc0c5865688dcf)
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 
791184866aSLeila Ghaffari // Wind Options for Advection
801184866aSLeila Ghaffari typedef enum {
811184866aSLeila Ghaffari   ADVECTION_WIND_ROTATION = 0,
821184866aSLeila Ghaffari   ADVECTION_WIND_TRANSLATION = 1,
831184866aSLeila Ghaffari } WindType;
841184866aSLeila Ghaffari static const char *const WindTypes[] = {
851184866aSLeila Ghaffari   "rotation",
861184866aSLeila Ghaffari   "translation",
871184866aSLeila Ghaffari   "WindType", "ADVECTION_WIND_", NULL
881184866aSLeila Ghaffari };
891184866aSLeila Ghaffari 
90ccaff030SJeremy L Thompson typedef enum {
91ccaff030SJeremy L Thompson   STAB_NONE = 0,
92ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
93ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
94ccaff030SJeremy L Thompson } StabilizationType;
95ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
9684d34d69SLeila Ghaffari   "none",
97ccaff030SJeremy L Thompson   "SU",
98ccaff030SJeremy L Thompson   "SUPG",
99ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
100ccaff030SJeremy L Thompson };
101ccaff030SJeremy L Thompson 
10284d34d69SLeila Ghaffari // Test Options
10384d34d69SLeila Ghaffari typedef enum {
10484d34d69SLeila Ghaffari   TEST_NONE = 0,               // Non test mode
10584d34d69SLeila Ghaffari   TEST_EXPLICIT = 1,           // Explicit test
10684d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab
10784d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab
10884d34d69SLeila Ghaffari } testType;
10984d34d69SLeila Ghaffari static const char *const testTypes[] = {
11084d34d69SLeila Ghaffari   "none",
11184d34d69SLeila Ghaffari   "explicit",
11284d34d69SLeila Ghaffari   "implicit_stab_none",
11384d34d69SLeila Ghaffari   "implicit_stab_supg",
11484d34d69SLeila Ghaffari   "testType", "TEST_", NULL
11584d34d69SLeila Ghaffari };
11684d34d69SLeila Ghaffari 
11784d34d69SLeila Ghaffari // Tests specific data
11884d34d69SLeila Ghaffari typedef struct {
11984d34d69SLeila Ghaffari   PetscScalar testtol;
12084d34d69SLeila Ghaffari   const char *filepath;
12184d34d69SLeila Ghaffari } testData;
12284d34d69SLeila Ghaffari 
12384d34d69SLeila Ghaffari testData testOptions[] = {
12484d34d69SLeila Ghaffari   [TEST_NONE] = {
12584d34d69SLeila Ghaffari     .testtol = 0.,
12684d34d69SLeila Ghaffari     .filepath = NULL
12784d34d69SLeila Ghaffari   },
12884d34d69SLeila Ghaffari   [TEST_EXPLICIT] = {
12984d34d69SLeila Ghaffari     .testtol = 1E-5,
13084d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin"
13184d34d69SLeila Ghaffari   },
13284d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_NONE] = {
13384d34d69SLeila Ghaffari     .testtol = 5E-4,
13484d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin"
13584d34d69SLeila Ghaffari   },
13684d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_SUPG] = {
13784d34d69SLeila Ghaffari     .testtol = 5E-4,
13884d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin"
13984d34d69SLeila Ghaffari   }
14084d34d69SLeila Ghaffari };
14184d34d69SLeila Ghaffari 
142ccaff030SJeremy L Thompson // Problem specific data
143ccaff030SJeremy L Thompson typedef struct {
1448b982baeSLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
145ebb4b9bdSLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction,
146ebb4b9bdSLeila Ghaffari                     applySur;
147ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
148ccaff030SJeremy L Thompson                        PetscScalar[], void *);
1497659d40cSLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
1507659d40cSLeila Ghaffari         *applyVol_ifunction_loc, *applySur_loc;
151ccaff030SJeremy L Thompson   const bool non_zero_time;
152ccaff030SJeremy L Thompson } problemData;
153ccaff030SJeremy L Thompson 
154ccaff030SJeremy L Thompson problemData problemOptions[] = {
155ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
156ccaff030SJeremy L Thompson     .dim                       = 3,
157ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1588b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
159b0137797SLeila Ghaffari     .setupVol                  = Setup,
160b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
161356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
162356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
163ccaff030SJeremy L Thompson     .ics                       = ICsDC,
164ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
165c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
166c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
167c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
168c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
169ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
17084d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
171ccaff030SJeremy L Thompson   },
172ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
173ccaff030SJeremy L Thompson     .dim                       = 3,
174ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1758b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
176b0137797SLeila Ghaffari     .setupVol                  = Setup,
177b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
178356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
179356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
180ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
181ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
182c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
183c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
184c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
185c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
1867659d40cSLeila Ghaffari     .applySur                  = Advection_Sur,
1877659d40cSLeila Ghaffari     .applySur_loc              = Advection_Sur_loc,
188ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
18984d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
190ccaff030SJeremy L Thompson   },
191ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
192ccaff030SJeremy L Thompson     .dim                       = 2,
193ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
1948b982baeSLeila Ghaffari     .qdatasizeSur              = 3,
195c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
196c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
197b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
198b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
199ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
200ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
201c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
202c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
203c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
204c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
2057659d40cSLeila Ghaffari     .applySur                  = Advection2d_Sur,
2067659d40cSLeila Ghaffari     .applySur_loc              = Advection2d_Sur_loc,
207ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
20884d34d69SLeila Ghaffari     .non_zero_time             = PETSC_TRUE,
209ccaff030SJeremy L Thompson   },
210ccaff030SJeremy L Thompson };
211ccaff030SJeremy L Thompson 
212ccaff030SJeremy L Thompson // PETSc user data
213ccaff030SJeremy L Thompson typedef struct User_ *User;
214ccaff030SJeremy L Thompson typedef struct Units_ *Units;
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson struct User_ {
217ccaff030SJeremy L Thompson   MPI_Comm comm;
218ccaff030SJeremy L Thompson   PetscInt outputfreq;
219ccaff030SJeremy L Thompson   DM dm;
220ccaff030SJeremy L Thompson   DM dmviz;
221ccaff030SJeremy L Thompson   Mat interpviz;
222ccaff030SJeremy L Thompson   Ceed ceed;
223ccaff030SJeremy L Thompson   Units units;
224ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
2251e150236SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction;
226ccaff030SJeremy L Thompson   Vec M;
227ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
228ccaff030SJeremy L Thompson   PetscInt contsteps;
229ccaff030SJeremy L Thompson };
230ccaff030SJeremy L Thompson 
231ccaff030SJeremy L Thompson struct Units_ {
232ccaff030SJeremy L Thompson   // fundamental units
233ccaff030SJeremy L Thompson   PetscScalar meter;
234ccaff030SJeremy L Thompson   PetscScalar kilogram;
235ccaff030SJeremy L Thompson   PetscScalar second;
236ccaff030SJeremy L Thompson   PetscScalar Kelvin;
237ccaff030SJeremy L Thompson   // derived units
238ccaff030SJeremy L Thompson   PetscScalar Pascal;
239ccaff030SJeremy L Thompson   PetscScalar JperkgK;
240ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
241ccaff030SJeremy L Thompson   PetscScalar WpermK;
242ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
243ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
244ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
24516c0476cSLeila Ghaffari   PetscScalar Joule;
246ccaff030SJeremy L Thompson };
247ccaff030SJeremy L Thompson 
248ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
249ccaff030SJeremy L Thompson struct SimpleBC_ {
2507659d40cSLeila Ghaffari   PetscInt nwall, nslip[3];
2517659d40cSLeila Ghaffari   PetscInt walls[6], slips[3][6];
25284d34d69SLeila Ghaffari   PetscBool userbc;
253ccaff030SJeremy L Thompson };
254ccaff030SJeremy L Thompson 
255ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
256ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
257ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
258ccaff030SJeremy L Thompson }
259ccaff030SJeremy L Thompson 
260ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
261ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
26284d34d69SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
263ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
264ccaff030SJeremy L Thompson 
265ccaff030SJeremy L Thompson   PetscSection section;
2661184866aSLeila Ghaffari   PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth;
2670c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2680c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
26984d34d69SLeila Ghaffari   Vec Uloc;
2700c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
271ccaff030SJeremy L Thompson   PetscErrorCode ierr;
272ccaff030SJeremy L Thompson 
273ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
274ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
275da51bdd9SLeila Ghaffari   dim -= height;
276ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
277ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
278ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
279ccaff030SJeremy L Thompson   fieldoff[0] = 0;
280ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
281ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
282ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
283ccaff030SJeremy L Thompson   }
284ccaff030SJeremy L Thompson 
2850c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2860c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2870c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2880c6c0b13SLeila Ghaffari   if (domainLabel) {
2890c6c0b13SLeila Ghaffari     IS domainIS;
2900c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
291*1119eeeeSJed Brown     if (domainIS) { // domainIS is non-empty
2920c6c0b13SLeila Ghaffari       ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2930c6c0b13SLeila Ghaffari       ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
294*1119eeeeSJed Brown     } else { // domainIS is NULL (empty)
295*1119eeeeSJed Brown       iterIS = NULL;
296*1119eeeeSJed Brown     }
2970c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2980c6c0b13SLeila Ghaffari   } else {
2990c6c0b13SLeila Ghaffari     iterIS = depthIS;
3000c6c0b13SLeila Ghaffari   }
301*1119eeeeSJed Brown   if (iterIS) {
3020c6c0b13SLeila Ghaffari     ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
3030c6c0b13SLeila Ghaffari     ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
304*1119eeeeSJed Brown   } else {
305*1119eeeeSJed Brown     Nelem = 0;
306*1119eeeeSJed Brown     iterIndices = NULL;
307*1119eeeeSJed Brown   }
308ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
3090c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
3100c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
311ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
31284d34d69SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
31384d34d69SLeila Ghaffari                                    &numindices, &indices, NULL, NULL);
31484d34d69SLeila Ghaffari     CHKERRQ(ierr);
31532b5ec5fSJed Brown     bool flip = false;
31632b5ec5fSJed Brown     if (height > 0) {
31732b5ec5fSJed Brown       PetscInt numCells, numFaces, start = -1;
31832b5ec5fSJed Brown       const PetscInt *orients, *faces, *cells;
31932b5ec5fSJed Brown       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
32032b5ec5fSJed Brown       ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr);
321ebb4b9bdSLeila Ghaffari       if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
322ebb4b9bdSLeila Ghaffari                                     "Expected one cell in support of exterior face, but got %D cells", numCells);
32332b5ec5fSJed Brown       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
32432b5ec5fSJed Brown       ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr);
32532b5ec5fSJed Brown       for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;}
326ebb4b9bdSLeila Ghaffari       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
327ebb4b9bdSLeila Ghaffari                                 "Could not find face %D in cone of its support", c);
32832b5ec5fSJed Brown       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
32932b5ec5fSJed Brown       if (orients[start] < 0) flip = true;
33032b5ec5fSJed Brown     }
33184d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
33284d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
33384d34d69SLeila Ghaffari           c);
334ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
335ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
33632b5ec5fSJed Brown       PetscInt ii = i;
33732b5ec5fSJed Brown       if (flip) {
33832b5ec5fSJed Brown         if (P == nnodes) ii = nnodes - 1 - i;
33932b5ec5fSJed Brown         else if (P*P == nnodes) {
34032b5ec5fSJed Brown           PetscInt row = i / P, col = i % P;
34132b5ec5fSJed Brown           ii = row + col * P;
342ebb4b9bdSLeila Ghaffari         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
343ebb4b9bdSLeila Ghaffari                           "No support for flipping point with %D nodes != P (%D) or P^2", nnodes, P);
34432b5ec5fSJed Brown       }
345ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
346ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
347ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
348ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
34932b5ec5fSJed Brown           if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j])
35032b5ec5fSJed Brown               != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j)
351ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
352ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
35332b5ec5fSJed Brown                      c, ii, f, j);
354ccaff030SJeremy L Thompson         }
355ccaff030SJeremy L Thompson       }
356ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
35732b5ec5fSJed Brown       PetscInt loc = Involute(indices[ii*ncomp[0]]);
3586f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
359ccaff030SJeremy L Thompson     }
36084d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
36184d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
36284d34d69SLeila Ghaffari     CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   }
3640c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3650c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3660c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
367ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
368*1119eeeeSJed Brown   if (iterIS) {
3690c6c0b13SLeila Ghaffari     ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
370*1119eeeeSJed Brown   }
3710c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3720c6c0b13SLeila Ghaffari 
373ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
374ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
375ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3766f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3776f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3786f55dfd5Svaleriabarra                             Erestrict);
379ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
380ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
381ccaff030SJeremy L Thompson }
382ccaff030SJeremy L Thompson 
383c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3841e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
3851e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
3861e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
3871e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
388c96c872fSLeila Ghaffari 
389c96c872fSLeila Ghaffari   DM dmcoord;
3901e150236SLeila Ghaffari   CeedInt dim, localNelem;
3911e150236SLeila Ghaffari   CeedInt Qdim;
392c96c872fSLeila Ghaffari   PetscErrorCode ierr;
393c96c872fSLeila Ghaffari 
394c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
3951e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
3961e150236SLeila Ghaffari   dim -= height;
3971e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
398c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
399c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
400c96c872fSLeila Ghaffari   CHKERRQ(ierr);
401ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value,
402ebb4b9bdSLeila Ghaffari                                    restrictq);
403c96c872fSLeila Ghaffari   CHKERRQ(ierr);
404ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value,
405ebb4b9bdSLeila Ghaffari                                    restrictx);
406c96c872fSLeila Ghaffari   CHKERRQ(ierr);
407c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
408c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
409c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
410c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
411c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
412c96c872fSLeila Ghaffari }
413c96c872fSLeila Ghaffari 
4141e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
4157659d40cSLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
4167659d40cSLeila Ghaffari     WindType wind_type, CeedOperator op_applyVol, CeedQFunction qf_applySur,
4177659d40cSLeila Ghaffari     CeedQFunction qf_setupSur, CeedInt height, CeedInt numP_Sur, CeedInt numQ_Sur,
4187659d40cSLeila Ghaffari     CeedInt qdatasizeSur, CeedInt NqptsSur, CeedBasis basisxSur,
4197659d40cSLeila Ghaffari     CeedBasis basisqSur, CeedOperator *op_apply) {
420ca3ac6ddSLeila Ghaffari 
4217659d40cSLeila Ghaffari   CeedInt dim, nFace;
4227659d40cSLeila Ghaffari   PetscInt lsize, localNelemSur[6];
4231e150236SLeila Ghaffari   Vec Xloc;
4247659d40cSLeila Ghaffari   CeedVector xcorners, qdataSur[6];
4257659d40cSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
4267659d40cSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
427ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
4281e150236SLeila Ghaffari   PetscScalar *x;
429ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
430ca3ac6ddSLeila Ghaffari 
431ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
432ca3ac6ddSLeila Ghaffari   // Composite Operaters
433ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
434ebb4b9bdSLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply,
435ebb4b9bdSLeila Ghaffari                               op_applyVol); // Apply a Sub-Operator for the volume
436ca3ac6ddSLeila Ghaffari 
4377659d40cSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION) {
4387659d40cSLeila Ghaffari     bc->nwall = 0;
4397659d40cSLeila Ghaffari     bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0;
4401e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
4411e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
4421e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
4431e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
4441e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
445ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
4467659d40cSLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
4477659d40cSLeila Ghaffari     if (dim == 2) nFace = 4;
4487659d40cSLeila Ghaffari     if (dim == 3) nFace = 6;
4499fe13df9SLeila Ghaffari 
4507659d40cSLeila Ghaffari     // Create CEED Operator for each boundary face
4517659d40cSLeila Ghaffari     for (CeedInt i=0; i<nFace; i++) {
4527659d40cSLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur,
4537659d40cSLeila Ghaffari                                      numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
4547659d40cSLeila Ghaffari                                      &restrictqdiSur[i]); CHKERRQ(ierr);
4557659d40cSLeila Ghaffari       // Create the CEED vectors that will be needed in Boundary setup
4567659d40cSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
4577659d40cSLeila Ghaffari       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
4587659d40cSLeila Ghaffari       // Create the operator that builds the quadrature data for the Boundary operator
4597659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
460ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur,
461ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4627659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
463ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
4647659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
465ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4667659d40cSLeila Ghaffari       // Create Boundary operator
4677659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
468ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur,
469ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4707659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
4717659d40cSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
4727659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, xcorners);
473ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur,
474ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4757659d40cSLeila Ghaffari       // Apply CEED operator for Boundary setup
476ebb4b9bdSLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i],
477ebb4b9bdSLeila Ghaffari                         CEED_REQUEST_IMMEDIATE);
4787659d40cSLeila Ghaffari       // Apply Sub-Operator for the Boundary
4797659d40cSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
4809fe13df9SLeila Ghaffari     }
4811e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
482ca3ac6ddSLeila Ghaffari   }
483ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
484ca3ac6ddSLeila Ghaffari }
485ca3ac6ddSLeila Ghaffari 
486ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
487ccaff030SJeremy L Thompson   PetscErrorCode ierr;
488ccaff030SJeremy L Thompson   PetscInt m;
489ccaff030SJeremy L Thompson 
490ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
491ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
494ccaff030SJeremy L Thompson }
495ccaff030SJeremy L Thompson 
496ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
497ccaff030SJeremy L Thompson   PetscErrorCode ierr;
498ccaff030SJeremy L Thompson   PetscInt mceed, mpetsc;
499ccaff030SJeremy L Thompson   PetscScalar *a;
500ccaff030SJeremy L Thompson 
501ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
502ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
503ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
504ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
50584d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
50684d34d69SLeila Ghaffari                                   mpetsc, mceed);
507ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
508ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
509ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
510ccaff030SJeremy L Thompson }
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
513ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
514ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
515ccaff030SJeremy L Thompson   PetscErrorCode ierr;
516ccaff030SJeremy L Thompson   Vec Qbc;
517ccaff030SJeremy L Thompson 
518ccaff030SJeremy L Thompson   PetscFunctionBegin;
519ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
520ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
521ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
522ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
523ccaff030SJeremy L Thompson }
524ccaff030SJeremy L Thompson 
525ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
526ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
527ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
528ccaff030SJeremy L Thompson   PetscErrorCode ierr;
529ccaff030SJeremy L Thompson   User user = *(User *)userData;
530ccaff030SJeremy L Thompson   PetscScalar *q, *g;
531ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
532ccaff030SJeremy L Thompson 
533ccaff030SJeremy L Thompson   // Global-to-local
534ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
535ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
536ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
538ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
539ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
540ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson 
543ccaff030SJeremy L Thompson   // Ceed Vectors
544ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
545ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
546ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
547ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
548ccaff030SJeremy L Thompson 
549ccaff030SJeremy L Thompson   // Apply CEED operator
550ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
551ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson   // Restore vectors
554ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
555ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
556ccaff030SJeremy L Thompson 
557ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson 
560ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
561ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
562ccaff030SJeremy L Thompson   CHKERRQ(ierr);
563ccaff030SJeremy L Thompson 
564ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
565ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
567ccaff030SJeremy L Thompson }
568ccaff030SJeremy L Thompson 
569ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
570ccaff030SJeremy L Thompson                                    void *userData) {
571ccaff030SJeremy L Thompson   PetscErrorCode ierr;
572ccaff030SJeremy L Thompson   User user = *(User *)userData;
573ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
574ccaff030SJeremy L Thompson   PetscScalar *g;
575ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
576ccaff030SJeremy L Thompson 
577ccaff030SJeremy L Thompson   // Global-to-local
578ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
579ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
580ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
581ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
582ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
584ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
585ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
586ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
589ccaff030SJeremy L Thompson 
590ccaff030SJeremy L Thompson   // Ceed Vectors
591ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
593ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
594ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
595ccaff030SJeremy L Thompson                      (PetscScalar *)q);
596ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
597ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
598ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
599ccaff030SJeremy L Thompson 
600ccaff030SJeremy L Thompson   // Apply CEED operator
601ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
602ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
603ccaff030SJeremy L Thompson 
604ccaff030SJeremy L Thompson   // Restore vectors
605ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
606ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
607ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson 
609ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
610ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
613ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
614ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
615ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
616ccaff030SJeremy L Thompson }
617ccaff030SJeremy L Thompson 
618ccaff030SJeremy L Thompson // User provided TS Monitor
619ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
620ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
621ccaff030SJeremy L Thompson   User user = ctx;
622ccaff030SJeremy L Thompson   Vec Qloc;
623ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
624ccaff030SJeremy L Thompson   PetscViewer viewer;
625ccaff030SJeremy L Thompson   PetscErrorCode ierr;
626ccaff030SJeremy L Thompson 
627ccaff030SJeremy L Thompson   // Set up output
628ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
629ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
630ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
631ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
632ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
633ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
634ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
635ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
636ccaff030SJeremy L Thompson 
637ccaff030SJeremy L Thompson   // Output
638ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
639ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
640ccaff030SJeremy L Thompson   CHKERRQ(ierr);
641ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
642ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
643ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
6449d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson   if (user->dmviz) {
646ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
647ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
648ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
649ccaff030SJeremy L Thompson 
650ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
651ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
652ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
653ccaff030SJeremy L Thompson     CHKERRQ(ierr);
654ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
656ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
657ccaff030SJeremy L Thompson     CHKERRQ(ierr);
658ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
659ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
660ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
661ccaff030SJeremy L Thompson     CHKERRQ(ierr);
662ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
663ccaff030SJeremy L Thompson                               filepath_refined,
664ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
665ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
666ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
667ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
668ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
669ccaff030SJeremy L Thompson   }
670ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
671ccaff030SJeremy L Thompson 
672ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
673ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
674ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
675ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
676ccaff030SJeremy L Thompson   CHKERRQ(ierr);
677ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
678ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
679ccaff030SJeremy L Thompson 
680ccaff030SJeremy L Thompson   // Save time stamp
681ccaff030SJeremy L Thompson   // Dimensionalize time back
682ccaff030SJeremy L Thompson   time /= user->units->second;
683ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
684ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
685ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
686ccaff030SJeremy L Thompson   CHKERRQ(ierr);
687ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
688ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
689ccaff030SJeremy L Thompson   #else
690ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
691ccaff030SJeremy L Thompson   #endif
692ccaff030SJeremy L Thompson   CHKERRQ(ierr);
693ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
694ccaff030SJeremy L Thompson 
695ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
696ccaff030SJeremy L Thompson }
697ccaff030SJeremy L Thompson 
69884d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
699ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
700ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
701ccaff030SJeremy L Thompson   PetscErrorCode ierr;
702ccaff030SJeremy L Thompson   CeedVector multlvec;
703ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
704ccaff030SJeremy L Thompson 
705ccaff030SJeremy L Thompson   ctxSetup->time = time;
706ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
707ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
708ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
709ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
710ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
711ccaff030SJeremy L Thompson 
712ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
713ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
714ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
715ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
716ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
717ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
718ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
719ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
720ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
721ccaff030SJeremy L Thompson   CHKERRQ(ierr);
722ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
723ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
724ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
725ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
726ccaff030SJeremy L Thompson 
727ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
728ccaff030SJeremy L Thompson }
729ccaff030SJeremy L Thompson 
730ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
731ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
732ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
733ccaff030SJeremy L Thompson   PetscErrorCode ierr;
734ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
735ccaff030SJeremy L Thompson   CeedOperator op_mass;
736ccaff030SJeremy L Thompson   CeedVector mceed;
737ccaff030SJeremy L Thompson   Vec Mloc;
738ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
739ccaff030SJeremy L Thompson 
740ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
741ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
742ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
743ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
744ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
745ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
746ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
747ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
748ccaff030SJeremy L Thompson 
749ccaff030SJeremy L Thompson   // Create the mass operator
750ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
751ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
752ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
753ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
754ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
755ccaff030SJeremy L Thompson 
756ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
757ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
758ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
759ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
760ccaff030SJeremy L Thompson 
761ccaff030SJeremy L Thompson   {
762ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
763ccaff030SJeremy L Thompson     CeedVector onesvec;
764ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
765ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
766ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
767ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
768ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
769ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
770ccaff030SJeremy L Thompson   }
771ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
772ccaff030SJeremy L Thompson 
773ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
774ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
775ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
776ccaff030SJeremy L Thompson 
777ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
778ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
779ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
780ccaff030SJeremy L Thompson }
781ccaff030SJeremy L Thompson 
78284d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
783ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
784ccaff030SJeremy L Thompson   PetscErrorCode ierr;
785ccaff030SJeremy L Thompson 
786ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
787ccaff030SJeremy L Thompson   {
788ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
789ccaff030SJeremy L Thompson     PetscFE fe;
790ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
791ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
792ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
79332ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
794ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
795ccaff030SJeremy L Thompson     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
796ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
79707af6069Svaleriabarra     {
79807af6069Svaleriabarra       PetscInt comps[1] = {1};
79907af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
80007af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
80107af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
80207af6069Svaleriabarra       comps[0] = 2;
80307af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
80407af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
80507af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
80607af6069Svaleriabarra       comps[0] = 3;
80707af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
80807af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
80907af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
81007af6069Svaleriabarra     }
81184d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
81284d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
81384d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
81484d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
81584d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
81684d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
81784d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
81884d34d69SLeila Ghaffari 
81984d34d69SLeila Ghaffari           }
82084d34d69SLeila Ghaffari         }
82184d34d69SLeila Ghaffari       }
82284d34d69SLeila Ghaffari     }
82384d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
82484d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
82584d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
82684d34d69SLeila Ghaffari     {
82784d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
82884d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
82984d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
83084d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
83184d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
83284d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
83384d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
83484d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
83584d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
83684d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
83784d34d69SLeila Ghaffari       } else
83884d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
83984d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
84084d34d69SLeila Ghaffari     }
841ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
842ccaff030SJeremy L Thompson     CHKERRQ(ierr);
843ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
844ccaff030SJeremy L Thompson   }
845ccaff030SJeremy L Thompson   {
846ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
847ccaff030SJeremy L Thompson     PetscSection section;
848ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
849ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
850ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
851ccaff030SJeremy L Thompson     CHKERRQ(ierr);
852ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
853ccaff030SJeremy L Thompson     CHKERRQ(ierr);
854ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
855ccaff030SJeremy L Thompson     CHKERRQ(ierr);
856ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
857ccaff030SJeremy L Thompson     CHKERRQ(ierr);
858ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
859ccaff030SJeremy L Thompson     CHKERRQ(ierr);
860ccaff030SJeremy L Thompson   }
861ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
862ccaff030SJeremy L Thompson }
863ccaff030SJeremy L Thompson 
864ccaff030SJeremy L Thompson int main(int argc, char **argv) {
865ccaff030SJeremy L Thompson   PetscInt ierr;
866ccaff030SJeremy L Thompson   MPI_Comm comm;
86784d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
868ccaff030SJeremy L Thompson   Mat interpviz;
869ccaff030SJeremy L Thompson   TS ts;
870ccaff030SJeremy L Thompson   TSAdapt adapt;
871ccaff030SJeremy L Thompson   User user;
872ccaff030SJeremy L Thompson   Units units;
873ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
87484d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
875ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
876ccaff030SJeremy L Thompson   PetscMPIInt rank;
877ccaff030SJeremy L Thompson   PetscScalar ftime;
878ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
879ccaff030SJeremy L Thompson   Ceed ceed;
88084d34d69SLeila Ghaffari   CeedInt numP, numQ;
881cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
88284d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
88384d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
884cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
885cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
886ccaff030SJeremy L Thompson   CeedScalar Rd;
88784d34d69SLeila Ghaffari   CeedMemType memtyperequested;
888ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
88916c0476cSLeila Ghaffari               kgpersquaredms, Joulepercubicm, Joule;
890ccaff030SJeremy L Thompson   problemType problemChoice;
891ccaff030SJeremy L Thompson   problemData *problem = NULL;
8921184866aSLeila Ghaffari   WindType wind_type;
893ccaff030SJeremy L Thompson   StabilizationType stab;
89484d34d69SLeila Ghaffari   testType testChoice;
89584d34d69SLeila Ghaffari   testData *test = NULL;
89684d34d69SLeila Ghaffari   PetscBool implicit;
897cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
898ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
89984d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
90084d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
901ccaff030SJeremy L Thompson   };
902ccaff030SJeremy L Thompson   double start, cpu_time_used;
90384d34d69SLeila Ghaffari   // Check PETSc CUDA support
90484d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
90584d34d69SLeila Ghaffari   // *INDENT-OFF*
90684d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
90784d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
90884d34d69SLeila Ghaffari   #else
90984d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
91084d34d69SLeila Ghaffari   #endif
91184d34d69SLeila Ghaffari   // *INDENT-ON*
912ccaff030SJeremy L Thompson 
913ccaff030SJeremy L Thompson   // Create the libCEED contexts
914ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
915ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
916ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
917ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
918ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
919ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
920ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
92116c0476cSLeila Ghaffari   CeedScalar E_wind      = 1.e6;     // J
922ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
923ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
924ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
925ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
926ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
927ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
928ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
929ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
930ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
931ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;       // [0,1]
932ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
933ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
934ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
935ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
936ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
937ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
938ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
939ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
940ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
94184d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
94284d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
943ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
9441184866aSLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0};
945ccaff030SJeremy L Thompson 
946ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
947ccaff030SJeremy L Thompson   if (ierr) return ierr;
948ccaff030SJeremy L Thompson 
949ccaff030SJeremy L Thompson   // Allocate PETSc context
950ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
951ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
952ccaff030SJeremy L Thompson 
953ccaff030SJeremy L Thompson   // Parse command line options
954ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
955ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
956ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
957ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
958ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
959ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
96084d34d69SLeila Ghaffari   testChoice = TEST_NONE;
96184d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
96284d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
96384d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
96484d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
96584d34d69SLeila Ghaffari   test = &testOptions[testChoice];
966ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
967ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
968ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
969ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
970ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
9711184866aSLeila Ghaffari   ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection",
9721184866aSLeila Ghaffari                           NULL, WindTypes, (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION),
9731184866aSLeila Ghaffari                           (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr);
9741184866aSLeila Ghaffari   PetscInt n = problem->dim;
97582c09b01SLeila Ghaffari   PetscBool userWind;
976ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsRealArray("-problem_advection_wind_translation",
977ebb4b9bdSLeila Ghaffari                                "Constant wind vector",
97882c09b01SLeila Ghaffari                                NULL, wind, &n, &userWind); CHKERRQ(ierr);
97982c09b01SLeila Ghaffari   if (wind_type == ADVECTION_WIND_ROTATION && userWind) {
980ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
981ebb4b9bdSLeila Ghaffari                        "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n");
98282c09b01SLeila Ghaffari     CHKERRQ(ierr);
9831184866aSLeila Ghaffari   }
984ebb4b9bdSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION
985ebb4b9bdSLeila Ghaffari       && problemChoice == NS_DENSITY_CURRENT) {
98667babd9cSLeila Ghaffari     SETERRQ(comm, PETSC_ERR_ARG_INCOMP,
98767babd9cSLeila Ghaffari             "-problem_advection_wind translation is not defined for -problem density_current");
98867babd9cSLeila Ghaffari   }
989ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
990ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
991ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
992ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
993ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
994ccaff030SJeremy L Thompson   CHKERRQ(ierr);
99584d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
99684d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
99784d34d69SLeila Ghaffari     CHKERRQ(ierr);
99884d34d69SLeila Ghaffari   }
999ccaff030SJeremy L Thompson   {
10007573aee6SJed Brown     PetscInt len;
10017573aee6SJed Brown     PetscBool flg;
1002ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
1003ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
1004ccaff030SJeremy L Thompson                                 NULL, bc.walls,
1005ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
1006ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
10077573aee6SJed Brown     if (flg) {
10087573aee6SJed Brown       bc.nwall = len;
10097573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
10107573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
10117573aee6SJed Brown     }
1012ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
1013ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1014ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
1015ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
1016ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
1017ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
1018ccaff030SJeremy L Thompson                                    &len), &flg);
1019ccaff030SJeremy L Thompson       CHKERRQ(ierr);
102084d34d69SLeila Ghaffari       if (flg) {
102184d34d69SLeila Ghaffari         bc.nslip[j] = len;
102284d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
102384d34d69SLeila Ghaffari       }
1024ccaff030SJeremy L Thompson     }
1025ccaff030SJeremy L Thompson   }
1026cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
1027cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
1028cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
1029ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1030ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1031ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1032ccaff030SJeremy L Thompson   meter = fabs(meter);
1033ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1034ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
1035ccaff030SJeremy L Thompson   second = fabs(second);
1036ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1037ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1038ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
1039ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
1040ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
1041ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1042ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
1043ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
1044ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
1045ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
1046ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
1047ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
1048ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
104916c0476cSLeila Ghaffari   ierr = PetscOptionsScalar("-E_wind", "Total energy of inflow wind",
105016c0476cSLeila Ghaffari                             NULL, E_wind, &E_wind, NULL); CHKERRQ(ierr);
1051ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
1052ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
1053ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1054ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1055ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1056ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1057ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1058ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1059ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1060ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1061ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1062ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1063ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1064ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1065ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1066ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1067ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1068ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
106984d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
107084d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
107184d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
107284d34d69SLeila Ghaffari     CHKERRQ(ierr);
107384d34d69SLeila Ghaffari   }
1074ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1075ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1076ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1077ccaff030SJeremy L Thompson   CHKERRQ(ierr);
107884d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
107984d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
108084d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
108184d34d69SLeila Ghaffari     CHKERRQ(ierr);
108284d34d69SLeila Ghaffari   }
1083ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1084ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1085ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1086ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1087ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1088ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1089ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1090ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1091ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1092ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1093ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1094ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1095ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1096ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
109782c09b01SLeila Ghaffari   n = problem->dim;
1098ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1099ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1100ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1101ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1102ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1103ccaff030SJeremy L Thompson   n = problem->dim;
1104ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1105ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1106ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1107ccaff030SJeremy L Thompson   {
1108ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1109ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1110ccaff030SJeremy L Thompson     if (norm > 0) {
1111ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1112ccaff030SJeremy L Thompson     }
1113ccaff030SJeremy L Thompson   }
1114ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1115ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1116ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1117ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1118ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
111984d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
112084d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
112184d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
112284d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
112381f92cf0SLeila Ghaffari   PetscBool userQextraSur;
1124ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsInt("-qextra_boundary",
1125ebb4b9bdSLeila Ghaffari                          "Number of extra quadrature points on in/outflow faces",
112681f92cf0SLeila Ghaffari                          NULL, qextraSur, &qextraSur, &userQextraSur); CHKERRQ(ierr);
1127ebb4b9bdSLeila Ghaffari   if ((wind_type == ADVECTION_WIND_ROTATION
1128ebb4b9bdSLeila Ghaffari        || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) {
1129ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
1130ebb4b9bdSLeila Ghaffari                        "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n");
113181f92cf0SLeila Ghaffari     CHKERRQ(ierr);
113281f92cf0SLeila Ghaffari   }
113384d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1134ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1135ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1136ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
113784d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
113884d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
113984d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
114084d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
114184d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
114284d34d69SLeila Ghaffari   CHKERRQ(ierr);
1143ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1144ccaff030SJeremy L Thompson 
1145ccaff030SJeremy L Thompson   // Define derived units
1146ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1147ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1148ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1149ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1150ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1151ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1152ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
115316c0476cSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
1154ccaff030SJeremy L Thompson 
1155ccaff030SJeremy L Thompson   // Scale variables to desired units
1156ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1157ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1158ccaff030SJeremy L Thompson   P0 *= Pascal;
115916c0476cSLeila Ghaffari   E_wind *= Joule;
1160ccaff030SJeremy L Thompson   N *= (1./second);
1161ccaff030SJeremy L Thompson   cv *= JperkgK;
1162ccaff030SJeremy L Thompson   cp *= JperkgK;
1163ccaff030SJeremy L Thompson   Rd = cp - cv;
1164ccaff030SJeremy L Thompson   g *= mpersquareds;
1165ccaff030SJeremy L Thompson   mu *= Pascal * second;
1166ccaff030SJeremy L Thompson   k *= WpermK;
1167ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1168ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1169ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1170ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1171ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1172ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1173ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1174ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1175ccaff030SJeremy L Thompson 
1176ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1177cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1178ccaff030SJeremy L Thompson   // Set up the libCEED context
1179ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1180ccaff030SJeremy L Thompson     .theta0 = theta0,
1181ccaff030SJeremy L Thompson     .thetaC = thetaC,
1182ccaff030SJeremy L Thompson     .P0 = P0,
1183ccaff030SJeremy L Thompson     .N = N,
1184ccaff030SJeremy L Thompson     .cv = cv,
1185ccaff030SJeremy L Thompson     .cp = cp,
1186ccaff030SJeremy L Thompson     .Rd = Rd,
1187ccaff030SJeremy L Thompson     .g = g,
1188ccaff030SJeremy L Thompson     .rc = rc,
1189ccaff030SJeremy L Thompson     .lx = lx,
1190ccaff030SJeremy L Thompson     .ly = ly,
1191ccaff030SJeremy L Thompson     .lz = lz,
1192ccaff030SJeremy L Thompson     .center[0] = center[0],
1193ccaff030SJeremy L Thompson     .center[1] = center[1],
1194ccaff030SJeremy L Thompson     .center[2] = center[2],
1195ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1196ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1197ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
11981184866aSLeila Ghaffari     .wind[0] = wind[0],
11991184866aSLeila Ghaffari     .wind[1] = wind[1],
12001184866aSLeila Ghaffari     .wind[2] = wind[2],
1201ccaff030SJeremy L Thompson     .time = 0,
12021184866aSLeila Ghaffari     .wind_type = wind_type,
1203ccaff030SJeremy L Thompson   };
1204ccaff030SJeremy L Thompson 
120584d34d69SLeila Ghaffari   // Create the mesh
1206ccaff030SJeremy L Thompson   {
1207ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1208ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
120984d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1210ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1211ccaff030SJeremy L Thompson   }
121284d34d69SLeila Ghaffari 
121384d34d69SLeila Ghaffari   // Distribute the mesh over processes
121484d34d69SLeila Ghaffari   {
1215ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1216ccaff030SJeremy L Thompson     PetscPartitioner part;
1217ccaff030SJeremy L Thompson 
1218ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1219ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1220ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1221ccaff030SJeremy L Thompson     if (dmDist) {
1222ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1223ccaff030SJeremy L Thompson       dm  = dmDist;
1224ccaff030SJeremy L Thompson     }
1225ccaff030SJeremy L Thompson   }
1226ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1227ccaff030SJeremy L Thompson 
122884d34d69SLeila Ghaffari   // Setup DM
1229ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1230ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
123184d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
123284d34d69SLeila Ghaffari 
123384d34d69SLeila Ghaffari   // Refine DM for high-order viz
1234ccaff030SJeremy L Thompson   dmviz = NULL;
1235ccaff030SJeremy L Thompson   interpviz = NULL;
1236ccaff030SJeremy L Thompson   if (viz_refine) {
1237ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1238ff6701fcSJed Brown 
1239ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1240ff6701fcSJed Brown     dmhierarchy[0] = dm;
124184d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1242ff6701fcSJed Brown       Mat interp_next;
1243ff6701fcSJed Brown 
1244ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1245ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1246bc6a41f7SJed Brown       ierr = DMClearDS(dmhierarchy[i+1]); CHKERRQ(ierr);
1247bc6a41f7SJed Brown       ierr = DMClearFields(dmhierarchy[i+1]); CHKERRQ(ierr);
1248ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1249ff6701fcSJed Brown       d = (d + 1) / 2;
1250ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1251ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1252ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1253ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1254ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1255ff6701fcSJed Brown       else {
1256ff6701fcSJed Brown         Mat C;
1257ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1258ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1259ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1260ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1261ff6701fcSJed Brown         interpviz = C;
1262ff6701fcSJed Brown       }
1263ff6701fcSJed Brown     }
1264cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1265ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1266cb3e2689Svaleriabarra     }
1267ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1268ccaff030SJeremy L Thompson   }
1269ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1270ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1271ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1272ccaff030SJeremy L Thompson   lnodes /= ncompq;
1273ccaff030SJeremy L Thompson 
127484d34d69SLeila Ghaffari   // Initialize CEED
127584d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
127684d34d69SLeila Ghaffari   // Set memtype
127784d34d69SLeila Ghaffari   CeedMemType memtypebackend;
127884d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
127984d34d69SLeila Ghaffari   // Check memtype compatibility
128084d34d69SLeila Ghaffari   if (!setmemtyperequest)
128184d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
128284d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
128384d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
128484d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
128584d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
128684d34d69SLeila Ghaffari 
128784d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
128884d34d69SLeila Ghaffari   numP = degree + 1;
128984d34d69SLeila Ghaffari   numQ = numP + qextra;
129084d34d69SLeila Ghaffari 
129184d34d69SLeila Ghaffari   // Print summary
129284d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1293ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1294ccaff030SJeremy L Thompson     int comm_size;
1295ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1296ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1297ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
129884d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1299ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1300ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1301ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
130284d34d69SLeila Ghaffari     const char *usedresource;
130384d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1304ccaff030SJeremy L Thompson 
130584d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
130684d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
130784d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
130884d34d69SLeila Ghaffari                        "  Problem:\n"
130984d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
131084d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
131184d34d69SLeila Ghaffari                        "  PETSc:\n"
131284d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
131384d34d69SLeila Ghaffari                        "  libCEED:\n"
131484d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
131584d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
131684d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
131784d34d69SLeila Ghaffari                        "  Mesh:\n"
131884d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
131984d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
132084d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
132184d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
132284d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
132384d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
132484d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
132584d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
132684d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
132784d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
132884d34d69SLeila Ghaffari                        (setmemtyperequest) ?
132984d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
133084d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
133184d34d69SLeila Ghaffari     CHKERRQ(ierr);
13320c6c0b13SLeila Ghaffari   }
13330c6c0b13SLeila Ghaffari 
1334ccaff030SJeremy L Thompson   // Set up global mass vector
1335ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1336ccaff030SJeremy L Thompson 
133784d34d69SLeila Ghaffari   // Set up libCEED
1338ccaff030SJeremy L Thompson   // CEED Bases
1339ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
134084d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
134184d34d69SLeila Ghaffari                                   &basisq);
134284d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
134384d34d69SLeila Ghaffari                                   &basisx);
134484d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
134584d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1346ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1347ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1348ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1349ccaff030SJeremy L Thompson 
1350ccaff030SJeremy L Thompson   // CEED Restrictions
13511e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
135284d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
135384d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1354ccaff030SJeremy L Thompson 
1355ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1356ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1357ccaff030SJeremy L Thompson 
1358ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1359bd910870SLeila Ghaffari   CeedInt NqptsVol;
136084d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
136184d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
13628b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
136384d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1364ccaff030SJeremy L Thompson 
1365ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1366ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1367ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1368ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1369ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
13708b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1371ccaff030SJeremy L Thompson 
1372ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1373ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1374ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1375ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1376ccaff030SJeremy L Thompson 
1377ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1378ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1379ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1380ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1381ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1382ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13838b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1384ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1385ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1386ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1387ccaff030SJeremy L Thompson   }
1388ccaff030SJeremy L Thompson 
1389ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1390ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1391ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1392ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1393ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1394ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1395ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
13968b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1397ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1398ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1399ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1400ccaff030SJeremy L Thompson   }
1401ccaff030SJeremy L Thompson 
1402ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1403ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
140484d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1405ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
140684d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
140784d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1408ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1409ccaff030SJeremy L Thompson 
1410ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1411ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
141284d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
141384d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1414ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1415ccaff030SJeremy L Thompson 
141684d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
141784d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
141884d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1419ccaff030SJeremy L Thompson 
1420ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1421ccaff030SJeremy L Thompson     CeedOperator op;
1422ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
142384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
142484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
142584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14268b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
142784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
142884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
142984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1430d3630711SJed Brown     user->op_rhs_vol = op;
1431ccaff030SJeremy L Thompson   }
1432ccaff030SJeremy L Thompson 
1433ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1434ccaff030SJeremy L Thompson     CeedOperator op;
1435ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
143684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
143784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
143884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
143984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14408b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
144184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
144284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
144384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1444d3630711SJed Brown     user->op_ifunction_vol = op;
1445ccaff030SJeremy L Thompson   }
1446ccaff030SJeremy L Thompson 
14476a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
14486a0edaf9SLeila Ghaffari   CeedInt height = 1;
14496a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
14501e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
14511e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1452cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1453cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1454cfa64770SLeila Ghaffari   CeedInt NqptsSur;
14557ed8b9c7SLeila Ghaffari   CeedQFunction qf_setupSur, qf_applySur;
1456cfa64770SLeila Ghaffari 
1457cfa64770SLeila Ghaffari   // CEED bases for the boundaries
1458ebb4b9bdSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur,
1459ebb4b9bdSLeila Ghaffari                                   CEED_GAUSS,
14606a0edaf9SLeila Ghaffari                                   &basisqSur);
14616a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
14626a0edaf9SLeila Ghaffari                                   &basisxSur);
14636a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
14646a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
14656a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
14666a0edaf9SLeila Ghaffari 
1467cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
14686a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
14696a0edaf9SLeila Ghaffari                               &qf_setupSur);
14706a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
14716a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
14726a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14736a0edaf9SLeila Ghaffari 
14747659d40cSLeila Ghaffari   // Creat Q-Function for Boundaries
14757ed8b9c7SLeila Ghaffari   qf_applySur = NULL;
14767659d40cSLeila Ghaffari   if (problem->applySur) {
14777659d40cSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur,
14787ed8b9c7SLeila Ghaffari                                 problem->applySur_loc, &qf_applySur);
14797ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "q", ncompq, CEED_EVAL_INTERP);
14807ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14817ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "x", ncompx, CEED_EVAL_INTERP);
14827ed8b9c7SLeila Ghaffari     CeedQFunctionAddOutput(qf_applySur, "v", ncompq, CEED_EVAL_INTERP);
14839fe13df9SLeila Ghaffari   }
14849fe13df9SLeila Ghaffari 
14859fe13df9SLeila Ghaffari   // Create CEED Operator for the whole domain
14869fe13df9SLeila Ghaffari   if (!implicit)
1487ebb4b9bdSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_rhs_vol,
1488ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
14897659d40cSLeila Ghaffari                                    height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
14907659d40cSLeila Ghaffari                                    basisqSur, &user->op_rhs); CHKERRQ(ierr);
14919fe13df9SLeila Ghaffari   if (implicit)
1492ebb4b9bdSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_ifunction_vol,
1493ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
14947659d40cSLeila Ghaffari                                    height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
14957659d40cSLeila Ghaffari                                    basisqSur, &user->op_ifunction); CHKERRQ(ierr);
14967659d40cSLeila Ghaffari   // Set up contex for QFunctions
1497ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1498ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1499ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1500ccaff030SJeremy L Thompson     .CtauS = CtauS,
1501ccaff030SJeremy L Thompson     .strong_form = strong_form,
1502ccaff030SJeremy L Thompson     .stabilization = stab,
1503ccaff030SJeremy L Thompson   };
15047659d40cSLeila Ghaffari   struct SurfaceContext_ ctxSurface = {
150516c0476cSLeila Ghaffari     .E_wind = E_wind,
15067659d40cSLeila Ghaffari     .strong_form = strong_form,
15077659d40cSLeila Ghaffari     .implicit = implicit,
15087659d40cSLeila Ghaffari   };
1509ccaff030SJeremy L Thompson   switch (problemChoice) {
1510ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1511ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1512ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1513ccaff030SJeremy L Thompson           sizeof ctxNS);
1514ccaff030SJeremy L Thompson     break;
1515ccaff030SJeremy L Thompson   case NS_ADVECTION:
1516ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1517ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1518ccaff030SJeremy L Thompson                                              sizeof ctxAdvection2d);
1519ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1520ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1521ebb4b9bdSLeila Ghaffari     if (qf_applySur) CeedQFunctionSetContext(qf_applySur, &ctxSurface,
1522ebb4b9bdSLeila Ghaffari           sizeof ctxSurface);
1523ccaff030SJeremy L Thompson   }
1524ccaff030SJeremy L Thompson 
1525ccaff030SJeremy L Thompson   // Set up PETSc context
1526ccaff030SJeremy L Thompson   // Set up units structure
1527ccaff030SJeremy L Thompson   units->meter = meter;
1528ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1529ccaff030SJeremy L Thompson   units->second = second;
1530ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1531ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1532ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1533ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1534ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1535ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1536ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1537ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
153816c0476cSLeila Ghaffari   units->Joule = Joule;
1539ccaff030SJeremy L Thompson 
1540ccaff030SJeremy L Thompson   // Set up user structure
1541ccaff030SJeremy L Thompson   user->comm = comm;
1542ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1543ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1544ccaff030SJeremy L Thompson   user->units = units;
1545ccaff030SJeremy L Thompson   user->dm = dm;
1546ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1547ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1548ccaff030SJeremy L Thompson   user->ceed = ceed;
1549ccaff030SJeremy L Thompson 
15508b982baeSLeila Ghaffari   // Calculate qdata and ICs
1551ccaff030SJeremy L Thompson   // Set up state global and local vectors
1552ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1553ccaff030SJeremy L Thompson 
1554cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1555ccaff030SJeremy L Thompson 
1556ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1557ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
15588b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
155984d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1560ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1561ccaff030SJeremy L Thompson 
156284d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
156384d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1564ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1565ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1566ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1567ccaff030SJeremy L Thompson     // slower execution.
1568ccaff030SJeremy L Thompson     Vec Qbc;
1569ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1570ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1571ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1572ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1573ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1574ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1575ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
157684d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
157784d34d69SLeila Ghaffari     CHKERRQ(ierr);
1578ccaff030SJeremy L Thompson   }
1579ccaff030SJeremy L Thompson 
1580ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1581ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1582ccaff030SJeremy L Thompson   // Gather initial Q values
1583ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1584ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1585ccaff030SJeremy L Thompson     PetscViewer viewer;
1586ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1587ccaff030SJeremy L Thompson     // Read input
1588ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1589ccaff030SJeremy L Thompson                          user->outputfolder);
1590ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1591ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1592ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1593ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1594ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1595ccaff030SJeremy L Thompson   }
1596ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1597ccaff030SJeremy L Thompson 
1598ccaff030SJeremy L Thompson // Create and setup TS
1599ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1600ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1601ccaff030SJeremy L Thompson   if (implicit) {
1602ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1603ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1604ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1605ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1606ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1607ccaff030SJeremy L Thompson     }
1608ccaff030SJeremy L Thompson   } else {
1609ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
1610ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1611ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1612ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1613ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1614ccaff030SJeremy L Thompson   }
1615ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1616ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1617ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
161884d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1619ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1620ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1621ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1622ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1623ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1624ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
162584d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1626ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1627ccaff030SJeremy L Thompson     }
1628ccaff030SJeremy L Thompson   } else { // continue from time of last output
1629ccaff030SJeremy L Thompson     PetscReal time;
1630ccaff030SJeremy L Thompson     PetscInt count;
1631ccaff030SJeremy L Thompson     PetscViewer viewer;
1632ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1633ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1634ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1635ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1636ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1637ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1638ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1639ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1640ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1641ccaff030SJeremy L Thompson   }
164284d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1643ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1644ccaff030SJeremy L Thompson   }
1645ccaff030SJeremy L Thompson 
1646ccaff030SJeremy L Thompson   // Solve
1647ccaff030SJeremy L Thompson   start = MPI_Wtime();
1648ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1649ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1650ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1651ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1652ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1653ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
165484d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1655ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
165684d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1657ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1658ccaff030SJeremy L Thompson   }
1659ccaff030SJeremy L Thompson 
1660ccaff030SJeremy L Thompson   // Get error
166184d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1662ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1663ccaff030SJeremy L Thompson     PetscReal norm;
1664ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1665ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1666ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1667ccaff030SJeremy L Thompson 
166884d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
166984d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1670ccaff030SJeremy L Thompson 
1671ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1672ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1673cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1674ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1675ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1676ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
167784d34d69SLeila Ghaffari     // Clean up vectors
167884d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
167984d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1680ccaff030SJeremy L Thompson   }
1681ccaff030SJeremy L Thompson 
1682ccaff030SJeremy L Thompson   // Output Statistics
1683ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr);
168484d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1685ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1686ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1687ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1688ccaff030SJeremy L Thompson   }
168984d34d69SLeila Ghaffari   // Output numerical values from command line
169084d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
169184d34d69SLeila Ghaffari 
169284d34d69SLeila Ghaffari   // compare reference solution values with current run
169384d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
169484d34d69SLeila Ghaffari     PetscViewer viewer;
169584d34d69SLeila Ghaffari     // Read reference file
169684d34d69SLeila Ghaffari     Vec Qref;
169784d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
169884d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
169984d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
170084d34d69SLeila Ghaffari     CHKERRQ(ierr);
170184d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
170284d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
170384d34d69SLeila Ghaffari 
170484d34d69SLeila Ghaffari     // Compute error with respect to reference solution
170584d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
170684d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
170784d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
170884d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
170984d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
171084d34d69SLeila Ghaffari     // Check error
171184d34d69SLeila Ghaffari     if (error > test->testtol) {
171284d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
171384d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
171484d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
171584d34d69SLeila Ghaffari     }
171684d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
171784d34d69SLeila Ghaffari   }
17189cf88b28Svaleriabarra 
1719ccaff030SJeremy L Thompson   // Clean up libCEED
17208b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1721ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1722ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1723ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1724ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
172584d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
172684d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
172784d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
172884d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
172984d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
173084d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1731ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1732ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1733ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1734ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1735ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1736ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
17376a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
17386a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
17396a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
17406a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
17416a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
17426a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
17436a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
17447ed8b9c7SLeila Ghaffari   CeedQFunctionDestroy(&qf_applySur);
1745ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1746ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1747ccaff030SJeremy L Thompson 
1748ccaff030SJeremy L Thompson   // Clean up PETSc
1749ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1750ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1751ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1752ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1753ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1754ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1755ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1756ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1757ccaff030SJeremy L Thompson   return PetscFinalize();
1758ccaff030SJeremy L Thompson }
1759ccaff030SJeremy L Thompson 
1760