xref: /libCEED/examples/fluids/navierstokes.c (revision f259b0549bd403572b5e0931b0b451ba8759bb85)
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);
2911119eeeeSJed Brown     if (domainIS) { // domainIS is non-empty
2920c6c0b13SLeila Ghaffari       ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2930c6c0b13SLeila Ghaffari       ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2941119eeeeSJed Brown     } else { // domainIS is NULL (empty)
2951119eeeeSJed Brown       iterIS = NULL;
2961119eeeeSJed Brown     }
2970c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2980c6c0b13SLeila Ghaffari   } else {
2990c6c0b13SLeila Ghaffari     iterIS = depthIS;
3000c6c0b13SLeila Ghaffari   }
3011119eeeeSJed Brown   if (iterIS) {
3020c6c0b13SLeila Ghaffari     ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
3030c6c0b13SLeila Ghaffari     ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3041119eeeeSJed Brown   } else {
3051119eeeeSJed Brown     Nelem = 0;
3061119eeeeSJed Brown     iterIndices = NULL;
3071119eeeeSJed 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,
322*f259b054Svaleriabarra                                     "Expected one cell in support of exterior face, but got %D cells",
323*f259b054Svaleriabarra                                     numCells);
32432b5ec5fSJed Brown       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
32532b5ec5fSJed Brown       ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr);
32632b5ec5fSJed Brown       for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;}
327ebb4b9bdSLeila Ghaffari       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
328*f259b054Svaleriabarra                                 "Could not find face %D in cone of its support",
329*f259b054Svaleriabarra                                 c);
33032b5ec5fSJed Brown       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
33132b5ec5fSJed Brown       if (orients[start] < 0) flip = true;
33232b5ec5fSJed Brown     }
33384d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
33484d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
33584d34d69SLeila Ghaffari           c);
336ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
337ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
33832b5ec5fSJed Brown       PetscInt ii = i;
33932b5ec5fSJed Brown       if (flip) {
34032b5ec5fSJed Brown         if (P == nnodes) ii = nnodes - 1 - i;
34132b5ec5fSJed Brown         else if (P*P == nnodes) {
34232b5ec5fSJed Brown           PetscInt row = i / P, col = i % P;
34332b5ec5fSJed Brown           ii = row + col * P;
344ebb4b9bdSLeila Ghaffari         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
345*f259b054Svaleriabarra                           "No support for flipping point with %D nodes != P (%D) or P^2",
346*f259b054Svaleriabarra                           nnodes, P);
34732b5ec5fSJed Brown       }
348ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
349ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
350ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
351ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
35232b5ec5fSJed Brown           if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j])
35332b5ec5fSJed Brown               != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j)
354ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
355ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
35632b5ec5fSJed Brown                      c, ii, f, j);
357ccaff030SJeremy L Thompson         }
358ccaff030SJeremy L Thompson       }
359ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
36032b5ec5fSJed Brown       PetscInt loc = Involute(indices[ii*ncomp[0]]);
3616f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
362ccaff030SJeremy L Thompson     }
36384d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
36484d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
36584d34d69SLeila Ghaffari     CHKERRQ(ierr);
366ccaff030SJeremy L Thompson   }
3670c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3680c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3690c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
370ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3711119eeeeSJed Brown   if (iterIS) {
3720c6c0b13SLeila Ghaffari     ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3731119eeeeSJed Brown   }
3740c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3750c6c0b13SLeila Ghaffari 
376ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
377ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
378ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3796f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3806f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3816f55dfd5Svaleriabarra                             Erestrict);
382ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
383ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
384ccaff030SJeremy L Thompson }
385ccaff030SJeremy L Thompson 
386c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3871e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
3881e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
3891e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
3901e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
391c96c872fSLeila Ghaffari 
392c96c872fSLeila Ghaffari   DM dmcoord;
3931e150236SLeila Ghaffari   CeedInt dim, localNelem;
3941e150236SLeila Ghaffari   CeedInt Qdim;
395c96c872fSLeila Ghaffari   PetscErrorCode ierr;
396c96c872fSLeila Ghaffari 
397c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
3981e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
3991e150236SLeila Ghaffari   dim -= height;
4001e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
401c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
402c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
403c96c872fSLeila Ghaffari   CHKERRQ(ierr);
404ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value,
405ebb4b9bdSLeila Ghaffari                                    restrictq);
406c96c872fSLeila Ghaffari   CHKERRQ(ierr);
407ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value,
408ebb4b9bdSLeila Ghaffari                                    restrictx);
409c96c872fSLeila Ghaffari   CHKERRQ(ierr);
410c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
411c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
412c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
413c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
414c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
415c96c872fSLeila Ghaffari }
416c96c872fSLeila Ghaffari 
4171e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
4187659d40cSLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
4197659d40cSLeila Ghaffari     WindType wind_type, CeedOperator op_applyVol, CeedQFunction qf_applySur,
4207659d40cSLeila Ghaffari     CeedQFunction qf_setupSur, CeedInt height, CeedInt numP_Sur, CeedInt numQ_Sur,
4217659d40cSLeila Ghaffari     CeedInt qdatasizeSur, CeedInt NqptsSur, CeedBasis basisxSur,
4227659d40cSLeila Ghaffari     CeedBasis basisqSur, CeedOperator *op_apply) {
423ca3ac6ddSLeila Ghaffari 
4247659d40cSLeila Ghaffari   CeedInt dim, nFace;
4257659d40cSLeila Ghaffari   PetscInt lsize, localNelemSur[6];
4261e150236SLeila Ghaffari   Vec Xloc;
4277659d40cSLeila Ghaffari   CeedVector xcorners, qdataSur[6];
4287659d40cSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
4297659d40cSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
430ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
4311e150236SLeila Ghaffari   PetscScalar *x;
432ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
433ca3ac6ddSLeila Ghaffari 
434ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
435ca3ac6ddSLeila Ghaffari   // Composite Operaters
436ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
437ebb4b9bdSLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply,
438ebb4b9bdSLeila Ghaffari                               op_applyVol); // Apply a Sub-Operator for the volume
439ca3ac6ddSLeila Ghaffari 
4407659d40cSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION) {
4417659d40cSLeila Ghaffari     bc->nwall = 0;
4427659d40cSLeila Ghaffari     bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0;
4431e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
4441e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
4451e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
4461e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
4471e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
448ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
4497659d40cSLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
4507659d40cSLeila Ghaffari     if (dim == 2) nFace = 4;
4517659d40cSLeila Ghaffari     if (dim == 3) nFace = 6;
4529fe13df9SLeila Ghaffari 
4537659d40cSLeila Ghaffari     // Create CEED Operator for each boundary face
4547659d40cSLeila Ghaffari     for (CeedInt i=0; i<nFace; i++) {
4557659d40cSLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur,
456*f259b054Svaleriabarra                                      numQ_Sur, qdatasizeSur, &restrictqSur[i],
457*f259b054Svaleriabarra                                      &restrictxSur[i], &restrictqdiSur[i]);
458*f259b054Svaleriabarra       CHKERRQ(ierr);
4597659d40cSLeila Ghaffari       // Create the CEED vectors that will be needed in Boundary setup
4607659d40cSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
461*f259b054Svaleriabarra       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur,
462*f259b054Svaleriabarra                        &qdataSur[i]);
4637659d40cSLeila Ghaffari       // Create the operator that builds the quadrature data for the Boundary operator
4647659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
465ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur,
466ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4677659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
468ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
4697659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
470ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4717659d40cSLeila Ghaffari       // Create Boundary operator
4727659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
473ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur,
474ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4757659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
4767659d40cSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
477*f259b054Svaleriabarra       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur,
478*f259b054Svaleriabarra                            xcorners);
479ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur,
480ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4817659d40cSLeila Ghaffari       // Apply CEED operator for Boundary setup
482ebb4b9bdSLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i],
483ebb4b9bdSLeila Ghaffari                         CEED_REQUEST_IMMEDIATE);
4847659d40cSLeila Ghaffari       // Apply Sub-Operator for the Boundary
4857659d40cSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
4869fe13df9SLeila Ghaffari     }
4871e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
488ca3ac6ddSLeila Ghaffari   }
489ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
490ca3ac6ddSLeila Ghaffari }
491ca3ac6ddSLeila Ghaffari 
492ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
493ccaff030SJeremy L Thompson   PetscErrorCode ierr;
494ccaff030SJeremy L Thompson   PetscInt m;
495ccaff030SJeremy L Thompson 
496ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
497ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
499ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
500ccaff030SJeremy L Thompson }
501ccaff030SJeremy L Thompson 
502ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
503ccaff030SJeremy L Thompson   PetscErrorCode ierr;
504ccaff030SJeremy L Thompson   PetscInt mceed, mpetsc;
505ccaff030SJeremy L Thompson   PetscScalar *a;
506ccaff030SJeremy L Thompson 
507ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
508ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
509ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
510ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
51184d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
51284d34d69SLeila Ghaffari                                   mpetsc, mceed);
513ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
514ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
515ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
516ccaff030SJeremy L Thompson }
517ccaff030SJeremy L Thompson 
518ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
519ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
520ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
521ccaff030SJeremy L Thompson   PetscErrorCode ierr;
522ccaff030SJeremy L Thompson   Vec Qbc;
523ccaff030SJeremy L Thompson 
524ccaff030SJeremy L Thompson   PetscFunctionBegin;
525ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
526ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
527ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
528ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
529ccaff030SJeremy L Thompson }
530ccaff030SJeremy L Thompson 
531ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
532ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
533ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
534ccaff030SJeremy L Thompson   PetscErrorCode ierr;
535ccaff030SJeremy L Thompson   User user = *(User *)userData;
536ccaff030SJeremy L Thompson   PetscScalar *q, *g;
537ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
538ccaff030SJeremy L Thompson 
539ccaff030SJeremy L Thompson   // Global-to-local
540ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
541ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
543ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
544ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
545ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
546ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
548ccaff030SJeremy L Thompson 
549ccaff030SJeremy L Thompson   // Ceed Vectors
550ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
553ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
554ccaff030SJeremy L Thompson 
555ccaff030SJeremy L Thompson   // Apply CEED operator
556ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
557ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
558ccaff030SJeremy L Thompson 
559ccaff030SJeremy L Thompson   // Restore vectors
560ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
562ccaff030SJeremy L Thompson 
563ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
565ccaff030SJeremy L Thompson 
566ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
567ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
568ccaff030SJeremy L Thompson   CHKERRQ(ierr);
569ccaff030SJeremy L Thompson 
570ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
571ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
572ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
573ccaff030SJeremy L Thompson }
574ccaff030SJeremy L Thompson 
575ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
576ccaff030SJeremy L Thompson                                    void *userData) {
577ccaff030SJeremy L Thompson   PetscErrorCode ierr;
578ccaff030SJeremy L Thompson   User user = *(User *)userData;
579ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
580ccaff030SJeremy L Thompson   PetscScalar *g;
581ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
582ccaff030SJeremy L Thompson 
583ccaff030SJeremy L Thompson   // Global-to-local
584ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
585ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
586ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
589ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
591ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
593ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
594ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
595ccaff030SJeremy L Thompson 
596ccaff030SJeremy L Thompson   // Ceed Vectors
597ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
598ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
599ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
600ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
601ccaff030SJeremy L Thompson                      (PetscScalar *)q);
602ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
603ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
604ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
605ccaff030SJeremy L Thompson 
606ccaff030SJeremy L Thompson   // Apply CEED operator
607ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
608ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
609ccaff030SJeremy L Thompson 
610ccaff030SJeremy L Thompson   // Restore vectors
611ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
612ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
613ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
614ccaff030SJeremy L Thompson 
615ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
616ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
617ccaff030SJeremy L Thompson 
618ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
619ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
620ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
621ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
622ccaff030SJeremy L Thompson }
623ccaff030SJeremy L Thompson 
624ccaff030SJeremy L Thompson // User provided TS Monitor
625ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
626ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
627ccaff030SJeremy L Thompson   User user = ctx;
628ccaff030SJeremy L Thompson   Vec Qloc;
629ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
630ccaff030SJeremy L Thompson   PetscViewer viewer;
631ccaff030SJeremy L Thompson   PetscErrorCode ierr;
632ccaff030SJeremy L Thompson 
633ccaff030SJeremy L Thompson   // Set up output
634ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
635ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
636ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
637ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
638ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
639ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
640ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
641ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson 
643ccaff030SJeremy L Thompson   // Output
644ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
645ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
646ccaff030SJeremy L Thompson   CHKERRQ(ierr);
647ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
648ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
649ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
6509d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
651ccaff030SJeremy L Thompson   if (user->dmviz) {
652ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
653ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
654ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
655ccaff030SJeremy L Thompson 
656ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
657ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
658ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
659ccaff030SJeremy L Thompson     CHKERRQ(ierr);
660ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
661ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
662ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
663ccaff030SJeremy L Thompson     CHKERRQ(ierr);
664ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
665ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
666ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
667ccaff030SJeremy L Thompson     CHKERRQ(ierr);
668ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
669ccaff030SJeremy L Thompson                               filepath_refined,
670ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
671ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
672ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
673ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
674ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
675ccaff030SJeremy L Thompson   }
676ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
677ccaff030SJeremy L Thompson 
678ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
679ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
680ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
681ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
682ccaff030SJeremy L Thompson   CHKERRQ(ierr);
683ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
684ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
685ccaff030SJeremy L Thompson 
686ccaff030SJeremy L Thompson   // Save time stamp
687ccaff030SJeremy L Thompson   // Dimensionalize time back
688ccaff030SJeremy L Thompson   time /= user->units->second;
689ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
690ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
691ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
692ccaff030SJeremy L Thompson   CHKERRQ(ierr);
693ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
694ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
695ccaff030SJeremy L Thompson   #else
696ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
697ccaff030SJeremy L Thompson   #endif
698ccaff030SJeremy L Thompson   CHKERRQ(ierr);
699ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
700ccaff030SJeremy L Thompson 
701ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
702ccaff030SJeremy L Thompson }
703ccaff030SJeremy L Thompson 
70484d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
705ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
706ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
707ccaff030SJeremy L Thompson   PetscErrorCode ierr;
708ccaff030SJeremy L Thompson   CeedVector multlvec;
709ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
710ccaff030SJeremy L Thompson 
711ccaff030SJeremy L Thompson   ctxSetup->time = time;
712ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
713ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
714ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
715ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
716ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
717ccaff030SJeremy L Thompson 
718ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
719ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
720ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
721ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
722ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
723ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
724ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
725ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
726ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
727ccaff030SJeremy L Thompson   CHKERRQ(ierr);
728ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
729ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
730ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
731ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
732ccaff030SJeremy L Thompson 
733ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
734ccaff030SJeremy L Thompson }
735ccaff030SJeremy L Thompson 
736ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
737ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
738ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
739ccaff030SJeremy L Thompson   PetscErrorCode ierr;
740ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
741ccaff030SJeremy L Thompson   CeedOperator op_mass;
742ccaff030SJeremy L Thompson   CeedVector mceed;
743ccaff030SJeremy L Thompson   Vec Mloc;
744ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
745ccaff030SJeremy L Thompson 
746ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
747ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
748ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
749ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
750ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
751ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
752ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
753ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
754ccaff030SJeremy L Thompson 
755ccaff030SJeremy L Thompson   // Create the mass operator
756ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
757ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
758ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
759ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
760ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
761ccaff030SJeremy L Thompson 
762ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
763ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
764ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
765ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson 
767ccaff030SJeremy L Thompson   {
768ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
769ccaff030SJeremy L Thompson     CeedVector onesvec;
770ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
771ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
772ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
773ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
774ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
775ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
776ccaff030SJeremy L Thompson   }
777ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
778ccaff030SJeremy L Thompson 
779ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
780ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
781ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
782ccaff030SJeremy L Thompson 
783ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
784ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
785ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
786ccaff030SJeremy L Thompson }
787ccaff030SJeremy L Thompson 
78884d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
789ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
790ccaff030SJeremy L Thompson   PetscErrorCode ierr;
791ccaff030SJeremy L Thompson 
792ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
793ccaff030SJeremy L Thompson   {
794ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
795ccaff030SJeremy L Thompson     PetscFE fe;
796ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
797ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
798ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
79932ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
800ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
801ccaff030SJeremy L Thompson     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
802ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
80307af6069Svaleriabarra     {
80407af6069Svaleriabarra       PetscInt comps[1] = {1};
80507af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
80607af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
80707af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
80807af6069Svaleriabarra       comps[0] = 2;
80907af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
81007af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
81107af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
81207af6069Svaleriabarra       comps[0] = 3;
81307af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
81407af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
81507af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
81607af6069Svaleriabarra     }
81784d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
81884d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
81984d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
82084d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
82184d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
82284d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
823*f259b054Svaleriabarra                        "Boundary condition already set on face %D!\n",
824*f259b054Svaleriabarra                        bc->walls[w]);
82584d34d69SLeila Ghaffari 
82684d34d69SLeila Ghaffari           }
82784d34d69SLeila Ghaffari         }
82884d34d69SLeila Ghaffari       }
82984d34d69SLeila Ghaffari     }
83084d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
83184d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
83284d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
83384d34d69SLeila Ghaffari     {
83484d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
83584d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
83684d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
83784d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
83884d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
83984d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
84084d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
84184d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
84284d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
84384d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
84484d34d69SLeila Ghaffari       } else
84584d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
84684d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
84784d34d69SLeila Ghaffari     }
848ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
849ccaff030SJeremy L Thompson     CHKERRQ(ierr);
850ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
851ccaff030SJeremy L Thompson   }
852ccaff030SJeremy L Thompson   {
853ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
854ccaff030SJeremy L Thompson     PetscSection section;
855ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
856ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
857ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
858ccaff030SJeremy L Thompson     CHKERRQ(ierr);
859ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
860ccaff030SJeremy L Thompson     CHKERRQ(ierr);
861ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
862ccaff030SJeremy L Thompson     CHKERRQ(ierr);
863ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
864ccaff030SJeremy L Thompson     CHKERRQ(ierr);
865ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
866ccaff030SJeremy L Thompson     CHKERRQ(ierr);
867ccaff030SJeremy L Thompson   }
868ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
869ccaff030SJeremy L Thompson }
870ccaff030SJeremy L Thompson 
871ccaff030SJeremy L Thompson int main(int argc, char **argv) {
872ccaff030SJeremy L Thompson   PetscInt ierr;
873ccaff030SJeremy L Thompson   MPI_Comm comm;
87484d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
875ccaff030SJeremy L Thompson   Mat interpviz;
876ccaff030SJeremy L Thompson   TS ts;
877ccaff030SJeremy L Thompson   TSAdapt adapt;
878ccaff030SJeremy L Thompson   User user;
879ccaff030SJeremy L Thompson   Units units;
880ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
88184d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
882ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
883ccaff030SJeremy L Thompson   PetscMPIInt rank;
884ccaff030SJeremy L Thompson   PetscScalar ftime;
885ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
886ccaff030SJeremy L Thompson   Ceed ceed;
88784d34d69SLeila Ghaffari   CeedInt numP, numQ;
888cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
88984d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
89084d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
891cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
892cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
893ccaff030SJeremy L Thompson   CeedScalar Rd;
89484d34d69SLeila Ghaffari   CeedMemType memtyperequested;
895ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
89616c0476cSLeila Ghaffari               kgpersquaredms, Joulepercubicm, Joule;
897ccaff030SJeremy L Thompson   problemType problemChoice;
898ccaff030SJeremy L Thompson   problemData *problem = NULL;
8991184866aSLeila Ghaffari   WindType wind_type;
900ccaff030SJeremy L Thompson   StabilizationType stab;
90184d34d69SLeila Ghaffari   testType testChoice;
90284d34d69SLeila Ghaffari   testData *test = NULL;
90384d34d69SLeila Ghaffari   PetscBool implicit;
904cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
905ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
90684d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
90784d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
908ccaff030SJeremy L Thompson   };
909ccaff030SJeremy L Thompson   double start, cpu_time_used;
91084d34d69SLeila Ghaffari   // Check PETSc CUDA support
91184d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
91284d34d69SLeila Ghaffari   // *INDENT-OFF*
91384d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
91484d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
91584d34d69SLeila Ghaffari   #else
91684d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
91784d34d69SLeila Ghaffari   #endif
91884d34d69SLeila Ghaffari   // *INDENT-ON*
919ccaff030SJeremy L Thompson 
920ccaff030SJeremy L Thompson   // Create the libCEED contexts
921ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
922ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
923ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
924ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
925ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
926ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
927ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
92816c0476cSLeila Ghaffari   CeedScalar E_wind      = 1.e6;     // J
929ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
930ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
931ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
932ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
933ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
934ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
935ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
936ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
937ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
938ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;       // [0,1]
939ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
940ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
941ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
942ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
943ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
944ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
945ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
946ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
947ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
94884d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
94984d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
950ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
9511184866aSLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0};
952ccaff030SJeremy L Thompson 
953ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
954ccaff030SJeremy L Thompson   if (ierr) return ierr;
955ccaff030SJeremy L Thompson 
956ccaff030SJeremy L Thompson   // Allocate PETSc context
957ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
958ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
959ccaff030SJeremy L Thompson 
960ccaff030SJeremy L Thompson   // Parse command line options
961ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
962ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
963ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
964ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
965ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
966ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
96784d34d69SLeila Ghaffari   testChoice = TEST_NONE;
96884d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
96984d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
97084d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
97184d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
97284d34d69SLeila Ghaffari   test = &testOptions[testChoice];
973ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
974ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
975ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
976ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
977ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
9781184866aSLeila Ghaffari   ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection",
979*f259b054Svaleriabarra                           NULL, WindTypes,
980*f259b054Svaleriabarra                           (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION),
9811184866aSLeila Ghaffari                           (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr);
9821184866aSLeila Ghaffari   PetscInt n = problem->dim;
98382c09b01SLeila Ghaffari   PetscBool userWind;
984ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsRealArray("-problem_advection_wind_translation",
985ebb4b9bdSLeila Ghaffari                                "Constant wind vector",
98682c09b01SLeila Ghaffari                                NULL, wind, &n, &userWind); CHKERRQ(ierr);
98782c09b01SLeila Ghaffari   if (wind_type == ADVECTION_WIND_ROTATION && userWind) {
988ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
989ebb4b9bdSLeila Ghaffari                        "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n");
99082c09b01SLeila Ghaffari     CHKERRQ(ierr);
9911184866aSLeila Ghaffari   }
992ebb4b9bdSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION
993ebb4b9bdSLeila Ghaffari       && problemChoice == NS_DENSITY_CURRENT) {
99467babd9cSLeila Ghaffari     SETERRQ(comm, PETSC_ERR_ARG_INCOMP,
99567babd9cSLeila Ghaffari             "-problem_advection_wind translation is not defined for -problem density_current");
99667babd9cSLeila Ghaffari   }
997ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
998ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
999ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
1000ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
1001ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
1002ccaff030SJeremy L Thompson   CHKERRQ(ierr);
100384d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
100484d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
100584d34d69SLeila Ghaffari     CHKERRQ(ierr);
100684d34d69SLeila Ghaffari   }
1007ccaff030SJeremy L Thompson   {
10087573aee6SJed Brown     PetscInt len;
10097573aee6SJed Brown     PetscBool flg;
1010ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
1011ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
1012ccaff030SJeremy L Thompson                                 NULL, bc.walls,
1013ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
1014ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
10157573aee6SJed Brown     if (flg) {
10167573aee6SJed Brown       bc.nwall = len;
10177573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
10187573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
10197573aee6SJed Brown     }
1020ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
1021ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1022ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
1023ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
1024ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
1025ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
1026ccaff030SJeremy L Thompson                                    &len), &flg);
1027ccaff030SJeremy L Thompson       CHKERRQ(ierr);
102884d34d69SLeila Ghaffari       if (flg) {
102984d34d69SLeila Ghaffari         bc.nslip[j] = len;
103084d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
103184d34d69SLeila Ghaffari       }
1032ccaff030SJeremy L Thompson     }
1033ccaff030SJeremy L Thompson   }
1034cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
1035cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
1036cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
1037ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1038ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1039ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1040ccaff030SJeremy L Thompson   meter = fabs(meter);
1041ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1042ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
1043ccaff030SJeremy L Thompson   second = fabs(second);
1044ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1045ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1046ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
1047ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
1048ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
1049ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1050ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
1051ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
1052ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
1053ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
1054ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
1055ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
1056ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
105716c0476cSLeila Ghaffari   ierr = PetscOptionsScalar("-E_wind", "Total energy of inflow wind",
105816c0476cSLeila Ghaffari                             NULL, E_wind, &E_wind, NULL); CHKERRQ(ierr);
1059ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
1060ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
1061ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1062ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1063ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1064ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1065ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1066ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1067ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1068ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1069ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1070ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1071ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1072ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1073ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1074ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1075ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1076ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
107784d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
107884d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
107984d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
108084d34d69SLeila Ghaffari     CHKERRQ(ierr);
108184d34d69SLeila Ghaffari   }
1082ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1083ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1084ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1085ccaff030SJeremy L Thompson   CHKERRQ(ierr);
108684d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
108784d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
108884d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
108984d34d69SLeila Ghaffari     CHKERRQ(ierr);
109084d34d69SLeila Ghaffari   }
1091ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1092ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1093ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1094ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1095ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1096ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1097ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1098ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1099ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1100ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1101ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1102ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1103ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1104ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
110582c09b01SLeila Ghaffari   n = problem->dim;
1106ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1107ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1108ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1109ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1110ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1111ccaff030SJeremy L Thompson   n = problem->dim;
1112ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1113ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1114ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1115ccaff030SJeremy L Thompson   {
1116ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1117ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1118ccaff030SJeremy L Thompson     if (norm > 0) {
1119ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1120ccaff030SJeremy L Thompson     }
1121ccaff030SJeremy L Thompson   }
1122ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1123ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1124ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1125ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1126ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
112784d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
112884d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
112984d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
113084d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
113181f92cf0SLeila Ghaffari   PetscBool userQextraSur;
1132ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsInt("-qextra_boundary",
1133ebb4b9bdSLeila Ghaffari                          "Number of extra quadrature points on in/outflow faces",
1134*f259b054Svaleriabarra                          NULL, qextraSur, &qextraSur, &userQextraSur);
1135*f259b054Svaleriabarra   CHKERRQ(ierr);
1136ebb4b9bdSLeila Ghaffari   if ((wind_type == ADVECTION_WIND_ROTATION
1137ebb4b9bdSLeila Ghaffari        || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) {
1138ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
1139ebb4b9bdSLeila Ghaffari                        "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n");
114081f92cf0SLeila Ghaffari     CHKERRQ(ierr);
114181f92cf0SLeila Ghaffari   }
114284d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1143ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1144ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1145ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
114684d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
114784d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
114884d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
114984d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
115084d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
115184d34d69SLeila Ghaffari   CHKERRQ(ierr);
1152ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1153ccaff030SJeremy L Thompson 
1154ccaff030SJeremy L Thompson   // Define derived units
1155ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1156ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1157ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1158ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1159ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1160ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1161ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
116216c0476cSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
1163ccaff030SJeremy L Thompson 
1164ccaff030SJeremy L Thompson   // Scale variables to desired units
1165ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1166ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1167ccaff030SJeremy L Thompson   P0 *= Pascal;
116816c0476cSLeila Ghaffari   E_wind *= Joule;
1169ccaff030SJeremy L Thompson   N *= (1./second);
1170ccaff030SJeremy L Thompson   cv *= JperkgK;
1171ccaff030SJeremy L Thompson   cp *= JperkgK;
1172ccaff030SJeremy L Thompson   Rd = cp - cv;
1173ccaff030SJeremy L Thompson   g *= mpersquareds;
1174ccaff030SJeremy L Thompson   mu *= Pascal * second;
1175ccaff030SJeremy L Thompson   k *= WpermK;
1176ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1177ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1178ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1179ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1180ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1181ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1182ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1183ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1184ccaff030SJeremy L Thompson 
1185ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1186cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1187ccaff030SJeremy L Thompson   // Set up the libCEED context
1188ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1189ccaff030SJeremy L Thompson     .theta0 = theta0,
1190ccaff030SJeremy L Thompson     .thetaC = thetaC,
1191ccaff030SJeremy L Thompson     .P0 = P0,
1192ccaff030SJeremy L Thompson     .N = N,
1193ccaff030SJeremy L Thompson     .cv = cv,
1194ccaff030SJeremy L Thompson     .cp = cp,
1195ccaff030SJeremy L Thompson     .Rd = Rd,
1196ccaff030SJeremy L Thompson     .g = g,
1197ccaff030SJeremy L Thompson     .rc = rc,
1198ccaff030SJeremy L Thompson     .lx = lx,
1199ccaff030SJeremy L Thompson     .ly = ly,
1200ccaff030SJeremy L Thompson     .lz = lz,
1201ccaff030SJeremy L Thompson     .center[0] = center[0],
1202ccaff030SJeremy L Thompson     .center[1] = center[1],
1203ccaff030SJeremy L Thompson     .center[2] = center[2],
1204ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1205ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1206ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
12071184866aSLeila Ghaffari     .wind[0] = wind[0],
12081184866aSLeila Ghaffari     .wind[1] = wind[1],
12091184866aSLeila Ghaffari     .wind[2] = wind[2],
1210ccaff030SJeremy L Thompson     .time = 0,
12111184866aSLeila Ghaffari     .wind_type = wind_type,
1212ccaff030SJeremy L Thompson   };
1213ccaff030SJeremy L Thompson 
121484d34d69SLeila Ghaffari   // Create the mesh
1215ccaff030SJeremy L Thompson   {
1216ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1217ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
121884d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1219ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1220ccaff030SJeremy L Thompson   }
122184d34d69SLeila Ghaffari 
122284d34d69SLeila Ghaffari   // Distribute the mesh over processes
122384d34d69SLeila Ghaffari   {
1224ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1225ccaff030SJeremy L Thompson     PetscPartitioner part;
1226ccaff030SJeremy L Thompson 
1227ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1228ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1229ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1230ccaff030SJeremy L Thompson     if (dmDist) {
1231ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1232ccaff030SJeremy L Thompson       dm  = dmDist;
1233ccaff030SJeremy L Thompson     }
1234ccaff030SJeremy L Thompson   }
1235ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1236ccaff030SJeremy L Thompson 
123784d34d69SLeila Ghaffari   // Setup DM
1238ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1239ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
124084d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
124184d34d69SLeila Ghaffari 
124284d34d69SLeila Ghaffari   // Refine DM for high-order viz
1243ccaff030SJeremy L Thompson   dmviz = NULL;
1244ccaff030SJeremy L Thompson   interpviz = NULL;
1245ccaff030SJeremy L Thompson   if (viz_refine) {
1246ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1247ff6701fcSJed Brown 
1248ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1249ff6701fcSJed Brown     dmhierarchy[0] = dm;
125084d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1251ff6701fcSJed Brown       Mat interp_next;
1252ff6701fcSJed Brown 
1253ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1254ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1255bc6a41f7SJed Brown       ierr = DMClearDS(dmhierarchy[i+1]); CHKERRQ(ierr);
1256bc6a41f7SJed Brown       ierr = DMClearFields(dmhierarchy[i+1]); CHKERRQ(ierr);
1257ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1258ff6701fcSJed Brown       d = (d + 1) / 2;
1259ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1260*f259b054Svaleriabarra       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup);
1261*f259b054Svaleriabarra       CHKERRQ(ierr);
1262ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1263ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1264ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1265ff6701fcSJed Brown       else {
1266ff6701fcSJed Brown         Mat C;
1267ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1268ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1269ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1270ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1271ff6701fcSJed Brown         interpviz = C;
1272ff6701fcSJed Brown       }
1273ff6701fcSJed Brown     }
1274cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1275ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1276cb3e2689Svaleriabarra     }
1277ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1278ccaff030SJeremy L Thompson   }
1279ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1280ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1281ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1282ccaff030SJeremy L Thompson   lnodes /= ncompq;
1283ccaff030SJeremy L Thompson 
128484d34d69SLeila Ghaffari   // Initialize CEED
128584d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
128684d34d69SLeila Ghaffari   // Set memtype
128784d34d69SLeila Ghaffari   CeedMemType memtypebackend;
128884d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
128984d34d69SLeila Ghaffari   // Check memtype compatibility
129084d34d69SLeila Ghaffari   if (!setmemtyperequest)
129184d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
129284d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
129384d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
129484d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
129584d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
129684d34d69SLeila Ghaffari 
129784d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
129884d34d69SLeila Ghaffari   numP = degree + 1;
129984d34d69SLeila Ghaffari   numQ = numP + qextra;
130084d34d69SLeila Ghaffari 
130184d34d69SLeila Ghaffari   // Print summary
130284d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1303ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1304ccaff030SJeremy L Thompson     int comm_size;
1305ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1306ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1307ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
130884d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1309ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1310ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1311ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
131284d34d69SLeila Ghaffari     const char *usedresource;
131384d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1314ccaff030SJeremy L Thompson 
131584d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
131684d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
131784d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
131884d34d69SLeila Ghaffari                        "  Problem:\n"
131984d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
132084d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
132184d34d69SLeila Ghaffari                        "  PETSc:\n"
132284d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
132384d34d69SLeila Ghaffari                        "  libCEED:\n"
132484d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
132584d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
132684d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
132784d34d69SLeila Ghaffari                        "  Mesh:\n"
132884d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
132984d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
133084d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
133184d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
133284d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
133384d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
133484d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
133584d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
133684d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
133784d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
133884d34d69SLeila Ghaffari                        (setmemtyperequest) ?
133984d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
134084d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
134184d34d69SLeila Ghaffari     CHKERRQ(ierr);
13420c6c0b13SLeila Ghaffari   }
13430c6c0b13SLeila Ghaffari 
1344ccaff030SJeremy L Thompson   // Set up global mass vector
1345ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1346ccaff030SJeremy L Thompson 
134784d34d69SLeila Ghaffari   // Set up libCEED
1348ccaff030SJeremy L Thompson   // CEED Bases
1349ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
135084d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
135184d34d69SLeila Ghaffari                                   &basisq);
135284d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
135384d34d69SLeila Ghaffari                                   &basisx);
135484d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
135584d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1356ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1357ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1358ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1359ccaff030SJeremy L Thompson 
1360ccaff030SJeremy L Thompson   // CEED Restrictions
13611e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
136284d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
136384d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1364ccaff030SJeremy L Thompson 
1365ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1366ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1367ccaff030SJeremy L Thompson 
1368ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1369bd910870SLeila Ghaffari   CeedInt NqptsVol;
137084d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
137184d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
13728b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
137384d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1374ccaff030SJeremy L Thompson 
1375ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1376ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1377ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1378ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1379ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
13808b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1381ccaff030SJeremy L Thompson 
1382ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1383ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1384ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1385ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1386ccaff030SJeremy L Thompson 
1387ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1388ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1389ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1390ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1391ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1392ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13938b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1394ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1395ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1396ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1397ccaff030SJeremy L Thompson   }
1398ccaff030SJeremy L Thompson 
1399ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1400ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1401ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1402ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1403ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1404ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1405ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
14068b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1407ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1408ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1409ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1410ccaff030SJeremy L Thompson   }
1411ccaff030SJeremy L Thompson 
1412ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1413ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
141484d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1415ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
141684d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
141784d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1418ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1419ccaff030SJeremy L Thompson 
1420ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1421ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
142284d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
142384d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1424ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1425ccaff030SJeremy L Thompson 
142684d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
142784d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
142884d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1429ccaff030SJeremy L Thompson 
1430ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1431ccaff030SJeremy L Thompson     CeedOperator op;
1432ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
143384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
143484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
143584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14368b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
143784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
143884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
143984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1440d3630711SJed Brown     user->op_rhs_vol = op;
1441ccaff030SJeremy L Thompson   }
1442ccaff030SJeremy L Thompson 
1443ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1444ccaff030SJeremy L Thompson     CeedOperator op;
1445ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
144684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
144784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
144884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
144984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14508b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
145184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
145284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
145384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1454d3630711SJed Brown     user->op_ifunction_vol = op;
1455ccaff030SJeremy L Thompson   }
1456ccaff030SJeremy L Thompson 
14576a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
14586a0edaf9SLeila Ghaffari   CeedInt height = 1;
14596a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
14601e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
14611e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1462cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1463cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1464cfa64770SLeila Ghaffari   CeedInt NqptsSur;
14657ed8b9c7SLeila Ghaffari   CeedQFunction qf_setupSur, qf_applySur;
1466cfa64770SLeila Ghaffari 
1467cfa64770SLeila Ghaffari   // CEED bases for the boundaries
1468ebb4b9bdSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur,
1469ebb4b9bdSLeila Ghaffari                                   CEED_GAUSS,
14706a0edaf9SLeila Ghaffari                                   &basisqSur);
14716a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
14726a0edaf9SLeila Ghaffari                                   &basisxSur);
14736a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
14746a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
14756a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
14766a0edaf9SLeila Ghaffari 
1477cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
14786a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
14796a0edaf9SLeila Ghaffari                               &qf_setupSur);
14806a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
14816a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
14826a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14836a0edaf9SLeila Ghaffari 
14847659d40cSLeila Ghaffari   // Creat Q-Function for Boundaries
14857ed8b9c7SLeila Ghaffari   qf_applySur = NULL;
14867659d40cSLeila Ghaffari   if (problem->applySur) {
14877659d40cSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur,
14887ed8b9c7SLeila Ghaffari                                 problem->applySur_loc, &qf_applySur);
14897ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "q", ncompq, CEED_EVAL_INTERP);
14907ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14917ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "x", ncompx, CEED_EVAL_INTERP);
14927ed8b9c7SLeila Ghaffari     CeedQFunctionAddOutput(qf_applySur, "v", ncompq, CEED_EVAL_INTERP);
14939fe13df9SLeila Ghaffari   }
14949fe13df9SLeila Ghaffari 
14959fe13df9SLeila Ghaffari   // Create CEED Operator for the whole domain
14969fe13df9SLeila Ghaffari   if (!implicit)
1497ebb4b9bdSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_rhs_vol,
1498ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
1499*f259b054Svaleriabarra                                    height, numP_Sur, numQ_Sur, qdatasizeSur,
1500*f259b054Svaleriabarra                                    NqptsSur, basisxSur, basisqSur,
1501*f259b054Svaleriabarra                                    &user->op_rhs); CHKERRQ(ierr);
15029fe13df9SLeila Ghaffari   if (implicit)
1503*f259b054Svaleriabarra     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type,
1504*f259b054Svaleriabarra                                    user->op_ifunction_vol,
1505ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
1506*f259b054Svaleriabarra                                    height, numP_Sur, numQ_Sur, qdatasizeSur,
1507*f259b054Svaleriabarra                                    NqptsSur, basisxSur, basisqSur,
1508*f259b054Svaleriabarra                                    &user->op_ifunction); CHKERRQ(ierr);
15097659d40cSLeila Ghaffari   // Set up contex for QFunctions
1510ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1511ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1512ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1513ccaff030SJeremy L Thompson     .CtauS = CtauS,
1514ccaff030SJeremy L Thompson     .strong_form = strong_form,
1515ccaff030SJeremy L Thompson     .stabilization = stab,
1516ccaff030SJeremy L Thompson   };
15177659d40cSLeila Ghaffari   struct SurfaceContext_ ctxSurface = {
151816c0476cSLeila Ghaffari     .E_wind = E_wind,
15197659d40cSLeila Ghaffari     .strong_form = strong_form,
15207659d40cSLeila Ghaffari     .implicit = implicit,
15217659d40cSLeila Ghaffari   };
1522ccaff030SJeremy L Thompson   switch (problemChoice) {
1523ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1524ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1525ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1526ccaff030SJeremy L Thompson           sizeof ctxNS);
1527ccaff030SJeremy L Thompson     break;
1528ccaff030SJeremy L Thompson   case NS_ADVECTION:
1529ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1530ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1531ccaff030SJeremy L Thompson                                              sizeof ctxAdvection2d);
1532ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1533ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1534ebb4b9bdSLeila Ghaffari     if (qf_applySur) CeedQFunctionSetContext(qf_applySur, &ctxSurface,
1535ebb4b9bdSLeila Ghaffari           sizeof ctxSurface);
1536ccaff030SJeremy L Thompson   }
1537ccaff030SJeremy L Thompson 
1538ccaff030SJeremy L Thompson   // Set up PETSc context
1539ccaff030SJeremy L Thompson   // Set up units structure
1540ccaff030SJeremy L Thompson   units->meter = meter;
1541ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1542ccaff030SJeremy L Thompson   units->second = second;
1543ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1544ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1545ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1546ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1547ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1548ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1549ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1550ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
155116c0476cSLeila Ghaffari   units->Joule = Joule;
1552ccaff030SJeremy L Thompson 
1553ccaff030SJeremy L Thompson   // Set up user structure
1554ccaff030SJeremy L Thompson   user->comm = comm;
1555ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1556ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1557ccaff030SJeremy L Thompson   user->units = units;
1558ccaff030SJeremy L Thompson   user->dm = dm;
1559ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1560ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1561ccaff030SJeremy L Thompson   user->ceed = ceed;
1562ccaff030SJeremy L Thompson 
15638b982baeSLeila Ghaffari   // Calculate qdata and ICs
1564ccaff030SJeremy L Thompson   // Set up state global and local vectors
1565ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1566ccaff030SJeremy L Thompson 
1567cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1568ccaff030SJeremy L Thompson 
1569ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1570ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
15718b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
157284d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1573ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1574ccaff030SJeremy L Thompson 
157584d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
157684d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1577ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1578ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1579ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1580ccaff030SJeremy L Thompson     // slower execution.
1581ccaff030SJeremy L Thompson     Vec Qbc;
1582ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1583ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1584ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1585ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1586ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1587ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1588ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
158984d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
159084d34d69SLeila Ghaffari     CHKERRQ(ierr);
1591ccaff030SJeremy L Thompson   }
1592ccaff030SJeremy L Thompson 
1593ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1594ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1595ccaff030SJeremy L Thompson   // Gather initial Q values
1596ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1597ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1598ccaff030SJeremy L Thompson     PetscViewer viewer;
1599ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1600ccaff030SJeremy L Thompson     // Read input
1601ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1602ccaff030SJeremy L Thompson                          user->outputfolder);
1603ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1604ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1605ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1607ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1608ccaff030SJeremy L Thompson   }
1609ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1610ccaff030SJeremy L Thompson 
1611ccaff030SJeremy L Thompson // Create and setup TS
1612ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1613ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1614ccaff030SJeremy L Thompson   if (implicit) {
1615ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1616ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1617ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1618ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1619ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1620ccaff030SJeremy L Thompson     }
1621ccaff030SJeremy L Thompson   } else {
1622ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
1623ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1624ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1625ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1626ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1627ccaff030SJeremy L Thompson   }
1628ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1629ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1630ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
163184d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1632ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1633ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1634ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1635ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1636ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1637ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
163884d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1639ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1640ccaff030SJeremy L Thompson     }
1641ccaff030SJeremy L Thompson   } else { // continue from time of last output
1642ccaff030SJeremy L Thompson     PetscReal time;
1643ccaff030SJeremy L Thompson     PetscInt count;
1644ccaff030SJeremy L Thompson     PetscViewer viewer;
1645ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1646ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1647ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1648ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1649ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1650ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1651ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1652ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1653ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1654ccaff030SJeremy L Thompson   }
165584d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1656ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1657ccaff030SJeremy L Thompson   }
1658ccaff030SJeremy L Thompson 
1659ccaff030SJeremy L Thompson   // Solve
1660ccaff030SJeremy L Thompson   start = MPI_Wtime();
1661ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1662ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1663ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1664ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1665ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1666ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
166784d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1668ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
166984d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1670ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1671ccaff030SJeremy L Thompson   }
1672ccaff030SJeremy L Thompson 
1673ccaff030SJeremy L Thompson   // Get error
167484d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1675ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1676ccaff030SJeremy L Thompson     PetscReal norm;
1677ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1678ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1679ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1680ccaff030SJeremy L Thompson 
168184d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
168284d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1683ccaff030SJeremy L Thompson 
1684ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1685ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1686cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1687ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1688ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1689ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
169084d34d69SLeila Ghaffari     // Clean up vectors
169184d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
169284d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1693ccaff030SJeremy L Thompson   }
1694ccaff030SJeremy L Thompson 
1695ccaff030SJeremy L Thompson   // Output Statistics
1696ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr);
169784d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1698ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1699ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1700ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1701ccaff030SJeremy L Thompson   }
170284d34d69SLeila Ghaffari   // Output numerical values from command line
170384d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
170484d34d69SLeila Ghaffari 
170584d34d69SLeila Ghaffari   // compare reference solution values with current run
170684d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
170784d34d69SLeila Ghaffari     PetscViewer viewer;
170884d34d69SLeila Ghaffari     // Read reference file
170984d34d69SLeila Ghaffari     Vec Qref;
171084d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
171184d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
171284d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
171384d34d69SLeila Ghaffari     CHKERRQ(ierr);
171484d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
171584d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
171684d34d69SLeila Ghaffari 
171784d34d69SLeila Ghaffari     // Compute error with respect to reference solution
171884d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
171984d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
172084d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
172184d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
172284d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
172384d34d69SLeila Ghaffari     // Check error
172484d34d69SLeila Ghaffari     if (error > test->testtol) {
172584d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
172684d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
172784d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
172884d34d69SLeila Ghaffari     }
172984d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
173084d34d69SLeila Ghaffari   }
17319cf88b28Svaleriabarra 
1732ccaff030SJeremy L Thompson   // Clean up libCEED
17338b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1734ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1735ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1736ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1737ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
173884d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
173984d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
174084d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
174184d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
174284d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
174384d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1744ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1745ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1746ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1747ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1748ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1749ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
17506a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
17516a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
17526a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
17536a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
17546a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
17556a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
17566a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
17577ed8b9c7SLeila Ghaffari   CeedQFunctionDestroy(&qf_applySur);
1758ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1759ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1760ccaff030SJeremy L Thompson 
1761ccaff030SJeremy L Thompson   // Clean up PETSc
1762ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1763ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1764ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1765ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1766ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1767ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1768ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1769ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1770ccaff030SJeremy L Thompson   return PetscFinalize();
1771ccaff030SJeremy L Thompson }
1772ccaff030SJeremy L Thompson 
1773