xref: /libCEED/examples/fluids/navierstokes.c (revision cfa6477014741e30cf537c865b8e05e33c245a9e)
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 //
31ea6e0f84SLeila Ghaffari //     ./navierstokes -ceed /cpu/self -problem density_current -degree_volume 1
32ea6e0f84SLeila Ghaffari //     ./navierstokes -ceed /gpu/occa -problem advection -degree_volume 1
33ccaff030SJeremy L Thompson //
34ea6e0f84SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test -degree_volume 1
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson /// @file
37ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc
38ccaff030SJeremy L Thompson 
39ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson #include <petscts.h>
42ccaff030SJeremy L Thompson #include <petscdmplex.h>
43ccaff030SJeremy L Thompson #include <ceed.h>
44ccaff030SJeremy L Thompson #include <stdbool.h>
45ccaff030SJeremy L Thompson #include <petscsys.h>
46ccaff030SJeremy L Thompson #include "common.h"
47b0137797SLeila Ghaffari #include "setup-boundary.h"
48ccaff030SJeremy L Thompson #include "advection.h"
49ccaff030SJeremy L Thompson #include "advection2d.h"
50ccaff030SJeremy L Thompson #include "densitycurrent.h"
51ccaff030SJeremy L Thompson 
52ccaff030SJeremy L Thompson // Problem Options
53ccaff030SJeremy L Thompson typedef enum {
54ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
55ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
56ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
57ccaff030SJeremy L Thompson } problemType;
58ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
59ccaff030SJeremy L Thompson   "density_current",
60ccaff030SJeremy L Thompson   "advection",
61ccaff030SJeremy L Thompson   "advection2d",
620c6c0b13SLeila Ghaffari   "problemType","NS_",0
63ccaff030SJeremy L Thompson };
64ccaff030SJeremy L Thompson 
65ccaff030SJeremy L Thompson typedef enum {
66ccaff030SJeremy L Thompson   STAB_NONE = 0,
67ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
68ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
69ccaff030SJeremy L Thompson } StabilizationType;
70ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
710c6c0b13SLeila Ghaffari   "NONE",
72ccaff030SJeremy L Thompson   "SU",
73ccaff030SJeremy L Thompson   "SUPG",
74ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
75ccaff030SJeremy L Thompson };
76ccaff030SJeremy L Thompson 
77ccaff030SJeremy L Thompson // Problem specific data
78ccaff030SJeremy L Thompson typedef struct {
798b982baeSLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
80ea6e0f84SLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applySur_rhs,
81ea6e0f84SLeila Ghaffari                     applyVol_ifunction, applySur_ifunction;
82ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
83ccaff030SJeremy L Thompson                        PetscScalar[], void *);
84ea6e0f84SLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
85ea6e0f84SLeila Ghaffari              *applySur_rhs_loc, *applyVol_ifunction_loc, *applySur_ifunction_loc;
86ccaff030SJeremy L Thompson   const bool non_zero_time;
87ccaff030SJeremy L Thompson } problemData;
88ccaff030SJeremy L Thompson 
89ccaff030SJeremy L Thompson problemData problemOptions[] = {
90ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
91ccaff030SJeremy L Thompson     .dim                       = 3,
92ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
938b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
94b0137797SLeila Ghaffari     .setupVol                  = Setup,
95b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
96356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
97356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
98ccaff030SJeremy L Thompson     .ics                       = ICsDC,
99ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
100c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
101c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
1028b982baeSLeila Ghaffari     .applySur_rhs              = DC_Sur,
1038b982baeSLeila Ghaffari     .applySur_rhs_loc          = DC_Sur_loc,
104c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
105c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
106ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_DC_Sur,
107ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_DC_Sur_loc,
108ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
1090c6c0b13SLeila Ghaffari     .non_zero_time             = false,
110ccaff030SJeremy L Thompson   },
111ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
112ccaff030SJeremy L Thompson     .dim                       = 3,
113ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1148b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
115b0137797SLeila Ghaffari     .setupVol                  = Setup,
116b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
117356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
118356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
119ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
120ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
121c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
122c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
12329448992SLeila Ghaffari     .applySur_rhs              = Advection_Sur,
12429448992SLeila Ghaffari     .applySur_rhs_loc          = Advection_Sur_loc,
125c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
126c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
127ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection_Sur,
128ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection_Sur_loc,
129ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
1300c6c0b13SLeila Ghaffari     .non_zero_time             = true,
131ccaff030SJeremy L Thompson   },
132ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
133ccaff030SJeremy L Thompson     .dim                       = 2,
134ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
1358b982baeSLeila Ghaffari     .qdatasizeSur              = 3,
136c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
137c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
138b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
139b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
140ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
141ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
142c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
143c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
144b0137797SLeila Ghaffari     .applySur_rhs              = Advection2d_Sur,
145b0137797SLeila Ghaffari     .applySur_rhs_loc          = Advection2d_Sur_loc,
146c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
147c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
148ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection2d_Sur,
149ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection2d_Sur_loc,
150ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
1510c6c0b13SLeila Ghaffari     .non_zero_time             = true,
152ccaff030SJeremy L Thompson   },
153ccaff030SJeremy L Thompson };
154ccaff030SJeremy L Thompson 
155ccaff030SJeremy L Thompson // PETSc user data
156ccaff030SJeremy L Thompson typedef struct User_ *User;
157ccaff030SJeremy L Thompson typedef struct Units_ *Units;
158ccaff030SJeremy L Thompson 
159ccaff030SJeremy L Thompson struct User_ {
160ccaff030SJeremy L Thompson   MPI_Comm comm;
161ccaff030SJeremy L Thompson   PetscInt outputfreq;
162ccaff030SJeremy L Thompson   DM dm;
163ccaff030SJeremy L Thompson   DM dmviz;
164ccaff030SJeremy L Thompson   Mat interpviz;
165ccaff030SJeremy L Thompson   Ceed ceed;
166ccaff030SJeremy L Thompson   Units units;
167ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
168*cfa64770SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs_sur[6], op_rhs,
169*cfa64770SLeila Ghaffari                op_ifunction_vol, op_ifunction_sur[6], op_ifunction;
170ccaff030SJeremy L Thompson   Vec M;
171ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
172ccaff030SJeremy L Thompson   PetscInt contsteps;
173ccaff030SJeremy L Thompson };
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson struct Units_ {
176ccaff030SJeremy L Thompson   // fundamental units
177ccaff030SJeremy L Thompson   PetscScalar meter;
178ccaff030SJeremy L Thompson   PetscScalar kilogram;
179ccaff030SJeremy L Thompson   PetscScalar second;
180ccaff030SJeremy L Thompson   PetscScalar Kelvin;
181ccaff030SJeremy L Thompson   // derived units
182ccaff030SJeremy L Thompson   PetscScalar Pascal;
183ccaff030SJeremy L Thompson   PetscScalar JperkgK;
184ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
185ccaff030SJeremy L Thompson   PetscScalar WpermK;
186ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
187ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
188ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
189ccaff030SJeremy L Thompson };
190ccaff030SJeremy L Thompson 
191ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
192ccaff030SJeremy L Thompson struct SimpleBC_ {
193*cfa64770SLeila Ghaffari   PetscInt nwall, nslip[3], noutflow;
194*cfa64770SLeila Ghaffari   PetscInt walls[10], slips[3][10], outflow[6];
195ccaff030SJeremy L Thompson };
196ccaff030SJeremy L Thompson 
197ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
198ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
199ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
200ccaff030SJeremy L Thompson }
201ccaff030SJeremy L Thompson 
202ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
203ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
204*cfa64770SLeila Ghaffari     CeedInt height, DMLabel domainLabel, PetscInt value,
205ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
206ccaff030SJeremy L Thompson 
207ccaff030SJeremy L Thompson   PetscSection   section;
2080c6c0b13SLeila Ghaffari   PetscInt       p, Nelem, Ndof, *erestrict, eoffset, nfields, dim,
2090c6c0b13SLeila Ghaffari                  depth;
2100c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2110c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
2120c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
213ccaff030SJeremy L Thompson   PetscErrorCode ierr;
214ccaff030SJeremy L Thompson   Vec Uloc;
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
217ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
218d1d4a8c6SJed Brown   dim -= height;
219ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm,&section); CHKERRQ(ierr);
220ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
221ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
222ccaff030SJeremy L Thompson   fieldoff[0] = 0;
223ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
224ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
225ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
226ccaff030SJeremy L Thompson   }
227ccaff030SJeremy L Thompson 
2280c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2290c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2300c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2310c6c0b13SLeila Ghaffari   if (domainLabel) {
2320c6c0b13SLeila Ghaffari     IS domainIS;
2330c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2340c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2350c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2360c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2370c6c0b13SLeila Ghaffari   } else {
2380c6c0b13SLeila Ghaffari     iterIS = depthIS;
2390c6c0b13SLeila Ghaffari   }
2400c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2410c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
242ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2430c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2440c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
245ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
2460c6c0b13SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
2470c6c0b13SLeila Ghaffari                                    &indices, NULL); CHKERRQ(ierr);
2480c6c0b13SLeila Ghaffari     if (numindices % fieldoff[nfields])
2490c6c0b13SLeila Ghaffari       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
2500c6c0b13SLeila Ghaffari                "Number of closure indices not compatible with Cell %D", c);
251ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
252ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
253ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
254ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
255ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
256ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
257ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
258ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
259ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
260ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
261ccaff030SJeremy L Thompson                      c, i, f, j);
262ccaff030SJeremy L Thompson         }
263ccaff030SJeremy L Thompson       }
264ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
265ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
2666f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
267ccaff030SJeremy L Thompson     }
2680c6c0b13SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
2690c6c0b13SLeila Ghaffari                                        &indices, NULL); CHKERRQ(ierr);
270ccaff030SJeremy L Thompson   }
2710c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
2720c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
2730c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
274ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
2750c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
2760c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
2770c6c0b13SLeila Ghaffari 
278ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
279ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
280ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
2816f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
2826f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
2836f55dfd5Svaleriabarra                             Erestrict);
284ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
285ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
286ccaff030SJeremy L Thompson }
287ccaff030SJeremy L Thompson 
288c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
289bd910870SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt ncompx, CeedInt dim,
290*cfa64770SLeila Ghaffari     CeedInt height, DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q,
291c96c872fSLeila Ghaffari     CeedInt qdatasize, CeedElemRestriction *restrictq,
292c96c872fSLeila Ghaffari     CeedElemRestriction *restrictx, CeedElemRestriction *restrictqdi) {
293c96c872fSLeila Ghaffari 
294c96c872fSLeila Ghaffari   DM dmcoord;
295c96c872fSLeila Ghaffari   CeedInt localNelem;
296c96c872fSLeila Ghaffari   CeedInt Qdim = CeedIntPow(Q, dim);
297c96c872fSLeila Ghaffari   PetscErrorCode ierr;
298c96c872fSLeila Ghaffari 
299c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
300c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
301c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
302c96c872fSLeila Ghaffari   CHKERRQ(ierr);
303c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
304c96c872fSLeila Ghaffari   CHKERRQ(ierr);
305c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
306c96c872fSLeila Ghaffari   CHKERRQ(ierr);
307c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
308c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
309c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
310c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
311c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
312c96c872fSLeila Ghaffari }
313c96c872fSLeila Ghaffari 
314ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
315ccaff030SJeremy L Thompson   PetscErrorCode ierr;
316ccaff030SJeremy L Thompson   PetscInt m;
317ccaff030SJeremy L Thompson 
318ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
319ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
320ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
321ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
322ccaff030SJeremy L Thompson }
323ccaff030SJeremy L Thompson 
324ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
325ccaff030SJeremy L Thompson   PetscErrorCode ierr;
326ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
327ccaff030SJeremy L Thompson   PetscScalar *a;
328ccaff030SJeremy L Thompson 
329ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
330ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
331ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
332ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
3330c6c0b13SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",mpetsc,mceed);
334ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
335ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
336ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
337ccaff030SJeremy L Thompson }
338ccaff030SJeremy L Thompson 
339ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
340ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
341ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
342ccaff030SJeremy L Thompson   PetscErrorCode ierr;
343ccaff030SJeremy L Thompson   Vec Qbc;
344ccaff030SJeremy L Thompson 
345ccaff030SJeremy L Thompson   PetscFunctionBegin;
346ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
347ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
348ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
349ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
350ccaff030SJeremy L Thompson }
351ccaff030SJeremy L Thompson 
352ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
353ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
354ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
355ccaff030SJeremy L Thompson   PetscErrorCode ierr;
356ccaff030SJeremy L Thompson   User user = *(User *)userData;
357ccaff030SJeremy L Thompson   PetscScalar *q, *g;
358ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
359ccaff030SJeremy L Thompson 
360ccaff030SJeremy L Thompson   // Global-to-local
361ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
362ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
365ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
366ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
367ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
368ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
369ccaff030SJeremy L Thompson 
370ccaff030SJeremy L Thompson   // Ceed Vectors
371ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
373ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
374ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
375ccaff030SJeremy L Thompson 
376ccaff030SJeremy L Thompson   // Apply CEED operator
377ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
378ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
379ccaff030SJeremy L Thompson 
380ccaff030SJeremy L Thompson   // Restore vectors
381ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
382ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
383ccaff030SJeremy L Thompson 
384ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
385ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
386ccaff030SJeremy L Thompson 
387ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
388ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
389ccaff030SJeremy L Thompson   CHKERRQ(ierr);
390ccaff030SJeremy L Thompson 
391ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
392ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
393ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
394ccaff030SJeremy L Thompson }
395ccaff030SJeremy L Thompson 
396ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
397ccaff030SJeremy L Thompson                                    void *userData) {
398ccaff030SJeremy L Thompson   PetscErrorCode ierr;
399ccaff030SJeremy L Thompson   User user = *(User *)userData;
400ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
401ccaff030SJeremy L Thompson   PetscScalar *g;
402ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
403ccaff030SJeremy L Thompson 
404ccaff030SJeremy L Thompson   // Global-to-local
405ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
406ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
408ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
409ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
411ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
412ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
413ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
414ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
415ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
416ccaff030SJeremy L Thompson 
417ccaff030SJeremy L Thompson   // Ceed Vectors
418ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
419ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
420ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
421ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
422ccaff030SJeremy L Thompson                      (PetscScalar *)q);
423ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
424ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
425ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
426ccaff030SJeremy L Thompson 
427ccaff030SJeremy L Thompson   // Apply CEED operator
428ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
429ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
430ccaff030SJeremy L Thompson 
431ccaff030SJeremy L Thompson   // Restore vectors
432ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
433ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
435ccaff030SJeremy L Thompson 
436ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
437ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
438ccaff030SJeremy L Thompson 
439ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
440ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
442ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
443ccaff030SJeremy L Thompson }
444ccaff030SJeremy L Thompson 
445ccaff030SJeremy L Thompson // User provided TS Monitor
446ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
447ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
448ccaff030SJeremy L Thompson   User user = ctx;
449ccaff030SJeremy L Thompson   Vec Qloc;
450ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
451ccaff030SJeremy L Thompson   PetscViewer viewer;
452ccaff030SJeremy L Thompson   PetscErrorCode ierr;
453ccaff030SJeremy L Thompson 
454ccaff030SJeremy L Thompson   // Set up output
455ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
456ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
457ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
458ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
459ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
460ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
461ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
462ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
463ccaff030SJeremy L Thompson 
464ccaff030SJeremy L Thompson   // Output
465ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
466ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
467ccaff030SJeremy L Thompson   CHKERRQ(ierr);
468ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
469ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
470ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
471ccaff030SJeremy L Thompson   if (user->dmviz) {
472ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
473ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
474ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
475ccaff030SJeremy L Thompson 
476ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
479ccaff030SJeremy L Thompson     CHKERRQ(ierr);
480ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
482ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
483ccaff030SJeremy L Thompson     CHKERRQ(ierr);
484ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
485ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
486ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
487ccaff030SJeremy L Thompson     CHKERRQ(ierr);
488ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
489ccaff030SJeremy L Thompson                               filepath_refined,
490ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
491ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
494ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
495ccaff030SJeremy L Thompson   }
4960c6c0b13SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson 
499ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
500ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
501ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
502ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
503ccaff030SJeremy L Thompson   CHKERRQ(ierr);
504ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
505ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
506ccaff030SJeremy L Thompson 
507ccaff030SJeremy L Thompson   // Save time stamp
508ccaff030SJeremy L Thompson   // Dimensionalize time back
509ccaff030SJeremy L Thompson   time /= user->units->second;
510ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
511ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
512ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
513ccaff030SJeremy L Thompson   CHKERRQ(ierr);
514ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
515ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
516ccaff030SJeremy L Thompson   #else
517ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
518ccaff030SJeremy L Thompson   #endif
519ccaff030SJeremy L Thompson   CHKERRQ(ierr);
520ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
521ccaff030SJeremy L Thompson 
522ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
523ccaff030SJeremy L Thompson }
524ccaff030SJeremy L Thompson 
5250c6c0b13SLeila Ghaffari static PetscErrorCode ICs_PetscMultiplicity(CeedOperator op_ics,
526ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
527ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
528ccaff030SJeremy L Thompson   PetscErrorCode ierr;
529ccaff030SJeremy L Thompson   CeedVector multlvec;
530ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
531ccaff030SJeremy L Thompson 
532ccaff030SJeremy L Thompson   ctxSetup->time = time;
533ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
536ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
538ccaff030SJeremy L Thompson 
539ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
540ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
542ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
543ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
544ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
545ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
546ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
548ccaff030SJeremy L Thompson   CHKERRQ(ierr);
549ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
550ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
553ccaff030SJeremy L Thompson 
554ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
555ccaff030SJeremy L Thompson }
556ccaff030SJeremy L Thompson 
557ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
558ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
559ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
560ccaff030SJeremy L Thompson   PetscErrorCode ierr;
561ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
562ccaff030SJeremy L Thompson   CeedOperator op_mass;
563ccaff030SJeremy L Thompson   CeedVector mceed;
564ccaff030SJeremy L Thompson   Vec Mloc;
565ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
566ccaff030SJeremy L Thompson 
567ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
568ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
569ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
570ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
571ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
572ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
573ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
574ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
575ccaff030SJeremy L Thompson 
576ccaff030SJeremy L Thompson   // Create the mass operator
577ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
578ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
579ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
580ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
581ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
582ccaff030SJeremy L Thompson 
583ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
584ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
585ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
586ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson 
588ccaff030SJeremy L Thompson   {
589ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
590ccaff030SJeremy L Thompson     CeedVector onesvec;
591ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
592ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
593ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
594ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
595ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
596ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
597ccaff030SJeremy L Thompson   }
598ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
599ccaff030SJeremy L Thompson 
600ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
601ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
602ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
603ccaff030SJeremy L Thompson 
604ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
605ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
606ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
607ccaff030SJeremy L Thompson }
608ccaff030SJeremy L Thompson 
6090c6c0b13SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
610ff6701fcSJed Brown                        SimpleBC bc, void *ctxSetup) {
611ccaff030SJeremy L Thompson   PetscErrorCode ierr;
612ccaff030SJeremy L Thompson 
613ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
614ccaff030SJeremy L Thompson   {
615ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
616ccaff030SJeremy L Thompson     PetscFE fe;
6170c6c0b13SLeila Ghaffari     PetscSpace fespace;
618ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
619ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
620ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
62132ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
622ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
623ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
624ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
6250c6c0b13SLeila Ghaffari     /* Wall boundary conditions are zero velocity and zero flux for density and energy */
6260c6c0b13SLeila Ghaffari     {
6270c6c0b13SLeila Ghaffari       PetscInt comps[3] = {1, 2, 3};
6280c6c0b13SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
6290c6c0b13SLeila Ghaffari                            3, comps, (void(*)(void))problem->bc,
6300c6c0b13SLeila Ghaffari                            bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
6310c6c0b13SLeila Ghaffari     }
63207af6069Svaleriabarra     {
63307af6069Svaleriabarra       PetscInt comps[1] = {1};
63407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
63507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
63607af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
63707af6069Svaleriabarra       comps[0] = 2;
63807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
63907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
64007af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
64107af6069Svaleriabarra       comps[0] = 3;
64207af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
64307af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
64407af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
64507af6069Svaleriabarra     }
646ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE,NULL);
647ccaff030SJeremy L Thompson     CHKERRQ(ierr);
6480c6c0b13SLeila Ghaffari     ierr = PetscFEGetBasisSpace(fe, &fespace); CHKERRQ(ierr);
649ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
650ccaff030SJeremy L Thompson   }
651ccaff030SJeremy L Thompson   {
652ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
653ccaff030SJeremy L Thompson     PetscSection section;
654ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
656ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
657ccaff030SJeremy L Thompson     CHKERRQ(ierr);
658ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
659ccaff030SJeremy L Thompson     CHKERRQ(ierr);
660ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
661ccaff030SJeremy L Thompson     CHKERRQ(ierr);
662ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
663ccaff030SJeremy L Thompson     CHKERRQ(ierr);
664ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
665ccaff030SJeremy L Thompson     CHKERRQ(ierr);
666ccaff030SJeremy L Thompson   }
667ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
668ccaff030SJeremy L Thompson }
669ccaff030SJeremy L Thompson 
670ccaff030SJeremy L Thompson int main(int argc, char **argv) {
671ccaff030SJeremy L Thompson   PetscInt ierr;
672ccaff030SJeremy L Thompson   MPI_Comm comm;
6730c6c0b13SLeila Ghaffari   DM dm, dmcoord, dmviz, dmvizfine;
674ccaff030SJeremy L Thompson   Mat interpviz;
675ccaff030SJeremy L Thompson   TS ts;
676ccaff030SJeremy L Thompson   TSAdapt adapt;
677ccaff030SJeremy L Thompson   User user;
678ccaff030SJeremy L Thompson   Units units;
679ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
680*cfa64770SLeila Ghaffari   PetscInt localNelemVol, lnodes, steps;
681ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
682ccaff030SJeremy L Thompson   PetscMPIInt rank;
683ccaff030SJeremy L Thompson   PetscScalar ftime;
684ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
685ccaff030SJeremy L Thompson   Ceed ceed;
686*cfa64770SLeila Ghaffari   CeedInt numP_Vol, numQ_Vol;
687*cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
688*cfa64770SLeila Ghaffari   CeedBasis basisxVol, basisxcVol, basisqVol;
689*cfa64770SLeila Ghaffari   CeedElemRestriction restrictxVol, restrictqVol, restrictqdiVol;
690*cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
691*cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
692ccaff030SJeremy L Thompson   CeedScalar Rd;
693ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
694ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
695ccaff030SJeremy L Thompson   problemType problemChoice;
696ccaff030SJeremy L Thompson   problemData *problem = NULL;
697ccaff030SJeremy L Thompson   StabilizationType stab;
6980c6c0b13SLeila Ghaffari   PetscBool   test, implicit;
699cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
700ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
7010c6c0b13SLeila Ghaffari     .nwall = 6,
7020c6c0b13SLeila Ghaffari     .walls = {1,2,3,4,5,6},
703ccaff030SJeremy L Thompson   };
704ccaff030SJeremy L Thompson   double start, cpu_time_used;
705ccaff030SJeremy L Thompson 
706ccaff030SJeremy L Thompson   // Create the libCEED contexts
707ccaff030SJeremy L Thompson   PetscScalar meter     = 1e-2;     // 1 meter in scaled length units
708ccaff030SJeremy L Thompson   PetscScalar second    = 1e-2;     // 1 second in scaled time units
709ccaff030SJeremy L Thompson   PetscScalar kilogram  = 1e-6;     // 1 kilogram in scaled mass units
710ccaff030SJeremy L Thompson   PetscScalar Kelvin    = 1;        // 1 Kelvin in scaled temperature units
711ccaff030SJeremy L Thompson   CeedScalar theta0     = 300.;     // K
712ccaff030SJeremy L Thompson   CeedScalar thetaC     = -15.;     // K
713ccaff030SJeremy L Thompson   CeedScalar P0         = 1.e5;     // Pa
714ccaff030SJeremy L Thompson   CeedScalar N          = 0.01;     // 1/s
715ccaff030SJeremy L Thompson   CeedScalar cv         = 717.;     // J/(kg K)
716ccaff030SJeremy L Thompson   CeedScalar cp         = 1004.;    // J/(kg K)
717ccaff030SJeremy L Thompson   CeedScalar g          = 9.81;     // m/s^2
718ccaff030SJeremy L Thompson   CeedScalar lambda     = -2./3.;   // -
719ccaff030SJeremy L Thompson   CeedScalar mu         = 75.;      // Pa s, dynamic viscosity
720ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
721ccaff030SJeremy L Thompson   CeedScalar k          = 0.02638;  // W/(m K)
722ccaff030SJeremy L Thompson   CeedScalar CtauS      = 0.;       // dimensionless
723ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
724ccaff030SJeremy L Thompson   PetscScalar lx        = 8000.;    // m
725ccaff030SJeremy L Thompson   PetscScalar ly        = 8000.;    // m
726ccaff030SJeremy L Thompson   PetscScalar lz        = 4000.;    // m
727ccaff030SJeremy L Thompson   CeedScalar rc         = 1000.;    // m (Radius of bubble)
728ccaff030SJeremy L Thompson   PetscScalar resx      = 1000.;    // m (resolution in x)
729ccaff030SJeremy L Thompson   PetscScalar resy      = 1000.;    // m (resolution in y)
730ccaff030SJeremy L Thompson   PetscScalar resz      = 1000.;    // m (resolution in z)
731ccaff030SJeremy L Thompson   PetscInt outputfreq   = 10;       // -
732ccaff030SJeremy L Thompson   PetscInt contsteps    = 0;        // -
733ea6e0f84SLeila Ghaffari   PetscInt degreeVol    = 1;        // -
734ea6e0f84SLeila Ghaffari   PetscInt degreeSur    = 1;        // -
735ea6e0f84SLeila Ghaffari   PetscInt qextraVol    = 2;        // -
736ea6e0f84SLeila Ghaffari   PetscInt qextraSur    = 2;        // -
7370c6c0b13SLeila Ghaffari   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
7380c6c0b13SLeila Ghaffari                                   DM_BOUNDARY_NONE
7390c6c0b13SLeila Ghaffari                                  };
740ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
741ccaff030SJeremy L Thompson 
742ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
743ccaff030SJeremy L Thompson   if (ierr) return ierr;
744ccaff030SJeremy L Thompson 
745ccaff030SJeremy L Thompson   // Allocate PETSc context
746ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
747ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
748ccaff030SJeremy L Thompson 
749ccaff030SJeremy L Thompson   // Parse command line options
750ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
751ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
752ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
753ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
754ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
755ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
7560c6c0b13SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
7570c6c0b13SLeila Ghaffari                           NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr);
758ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
759ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
760ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
761ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
762ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
763ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
764ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
765ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
767ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
768ccaff030SJeremy L Thompson   CHKERRQ(ierr);
769ccaff030SJeremy L Thompson   {
770*cfa64770SLeila Ghaffari     PetscInt len, len1;
771*cfa64770SLeila Ghaffari     PetscBool flg, flg1;
772ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
773ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
774ccaff030SJeremy L Thompson                                 NULL, bc.walls,
775ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
776ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
777ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
778ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
779ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
780ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
781ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
782ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
783ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
784ccaff030SJeremy L Thompson                                    &len), &flg);
785ccaff030SJeremy L Thompson       CHKERRQ(ierr);
7860c6c0b13SLeila Ghaffari       if (flg) bc.nslip[j] = len;
787ccaff030SJeremy L Thompson     }
788*cfa64770SLeila Ghaffari     ierr = PetscOptionsIntArray("-bc_outflow",
789*cfa64770SLeila Ghaffari                               "Use outflow boundary conditions on this list of faces",
790*cfa64770SLeila Ghaffari                               NULL, bc.outflow,
791*cfa64770SLeila Ghaffari                               (len1 = sizeof(bc.outflow) / sizeof(bc.outflow[0]),
792*cfa64770SLeila Ghaffari                               &len1), &flg1); CHKERRQ(ierr);
793*cfa64770SLeila Ghaffari     if (flg1) bc.noutflow = len1;
794ccaff030SJeremy L Thompson   }
795cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
796cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
797cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
798ccaff030SJeremy L Thompson   CHKERRQ(ierr);
799ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
800ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
801ccaff030SJeremy L Thompson   meter = fabs(meter);
802ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
803ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
804ccaff030SJeremy L Thompson   second = fabs(second);
805ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
806ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
807ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
808ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
809ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
810ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
811ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
812ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
813ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
814ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
815ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
816ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
817ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
818ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
819ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
820ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
821ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
822ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
823ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
824ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
825ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
826ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
827ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
828ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
829ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
830ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
831ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
832ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
833ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
834ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
835ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
836ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
837ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
838ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
839ccaff030SJeremy L Thompson   CHKERRQ(ierr);
840ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
841ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
842ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
843ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
844ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
845ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
846ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
847ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
848ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
849ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
850ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
851ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
852ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
853ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
854ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
8550c6c0b13SLeila Ghaffari   ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction",
8560c6c0b13SLeila Ghaffari                                NULL, DMBoundaryTypes, (PetscEnum *)periodicity,
8570c6c0b13SLeila Ghaffari                                &n, NULL); CHKERRQ(ierr);
8580c6c0b13SLeila Ghaffari   n = problem->dim;
859ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
860ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
861ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
862ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
863ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
864ccaff030SJeremy L Thompson   n = problem->dim;
865ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
866ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
867ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
868ccaff030SJeremy L Thompson   {
869ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
870ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
871ccaff030SJeremy L Thompson     if (norm > 0) {
872ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
873ccaff030SJeremy L Thompson     }
874ccaff030SJeremy L Thompson   }
875ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
876ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
877ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
878ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
879ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
880ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-degree_volume", "Polynomial degree of finite elements for the Volume",
881ea6e0f84SLeila Ghaffari                          NULL, degreeVol, &degreeVol, NULL); CHKERRQ(ierr);
882ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-qextra_volume", "Number of extra quadrature points for the Volume",
883ea6e0f84SLeila Ghaffari                          NULL, qextraVol, &qextraVol, NULL); CHKERRQ(ierr);
8840c6c0b13SLeila Ghaffari   PetscStrncpy(user->outputfolder, ".", 2);
885ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
886ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
887ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
888ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
889ccaff030SJeremy L Thompson 
890ccaff030SJeremy L Thompson   // Define derived units
891ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
892ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
893ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
894ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
895ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
896ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
897ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
898ccaff030SJeremy L Thompson 
899ccaff030SJeremy L Thompson   // Scale variables to desired units
900ccaff030SJeremy L Thompson   theta0 *= Kelvin;
901ccaff030SJeremy L Thompson   thetaC *= Kelvin;
902ccaff030SJeremy L Thompson   P0 *= Pascal;
903ccaff030SJeremy L Thompson   N *= (1./second);
904ccaff030SJeremy L Thompson   cv *= JperkgK;
905ccaff030SJeremy L Thompson   cp *= JperkgK;
906ccaff030SJeremy L Thompson   Rd = cp - cv;
907ccaff030SJeremy L Thompson   g *= mpersquareds;
908ccaff030SJeremy L Thompson   mu *= Pascal * second;
909ccaff030SJeremy L Thompson   k *= WpermK;
910ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
911ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
912ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
913ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
914ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
915ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
916ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
917ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
918ccaff030SJeremy L Thompson 
919ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
920*cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
921ccaff030SJeremy L Thompson   // Set up the libCEED context
922ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup =  {
923ccaff030SJeremy L Thompson     .theta0 = theta0,
924ccaff030SJeremy L Thompson     .thetaC = thetaC,
925ccaff030SJeremy L Thompson     .P0 = P0,
926ccaff030SJeremy L Thompson     .N = N,
927ccaff030SJeremy L Thompson     .cv = cv,
928ccaff030SJeremy L Thompson     .cp = cp,
929ccaff030SJeremy L Thompson     .Rd = Rd,
930ccaff030SJeremy L Thompson     .g = g,
931ccaff030SJeremy L Thompson     .rc = rc,
932ccaff030SJeremy L Thompson     .lx = lx,
933ccaff030SJeremy L Thompson     .ly = ly,
934ccaff030SJeremy L Thompson     .lz = lz,
9350c6c0b13SLeila Ghaffari     .periodicity0 = periodicity[0],
9360c6c0b13SLeila Ghaffari     .periodicity1 = periodicity[1],
9370c6c0b13SLeila Ghaffari     .periodicity2 = periodicity[2],
938ccaff030SJeremy L Thompson     .center[0] = center[0],
939ccaff030SJeremy L Thompson     .center[1] = center[1],
940ccaff030SJeremy L Thompson     .center[2] = center[2],
941ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
942ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
943ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
944ccaff030SJeremy L Thompson     .time = 0,
945ccaff030SJeremy L Thompson   };
946ccaff030SJeremy L Thompson 
947ccaff030SJeremy L Thompson   {
948ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
949ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
9500c6c0b13SLeila Ghaffari                                periodicity, PETSC_TRUE, &dm);
951ccaff030SJeremy L Thompson     CHKERRQ(ierr);
952ccaff030SJeremy L Thompson   }
9530c6c0b13SLeila Ghaffari   if (1) {
954ccaff030SJeremy L Thompson     DM               dmDist = NULL;
955ccaff030SJeremy L Thompson     PetscPartitioner part;
956ccaff030SJeremy L Thompson 
957ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
958ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
959ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
960ccaff030SJeremy L Thompson     if (dmDist) {
961ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
962ccaff030SJeremy L Thompson       dm  = dmDist;
963ccaff030SJeremy L Thompson     }
964ccaff030SJeremy L Thompson   }
965ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
966ccaff030SJeremy L Thompson 
967ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
968ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
969ea6e0f84SLeila Ghaffari   ierr = SetUpDM(dm, problem, degreeVol, &bc, &ctxSetup); CHKERRQ(ierr);
9700c6c0b13SLeila Ghaffari   if (!test) {
9710c6c0b13SLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
972ea6e0f84SLeila Ghaffari                        "Degree of Volumetric FEM Space: %D\n",
973ea6e0f84SLeila Ghaffari                        degreeVol); CHKERRQ(ierr);
9740c6c0b13SLeila Ghaffari   }
975ccaff030SJeremy L Thompson   dmviz = NULL;
976ccaff030SJeremy L Thompson   interpviz = NULL;
977ccaff030SJeremy L Thompson   if (viz_refine) {
978ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
979ff6701fcSJed Brown 
980ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
981ff6701fcSJed Brown     dmhierarchy[0] = dm;
982ea6e0f84SLeila Ghaffari     for (PetscInt i = 0, d = degreeVol; i < viz_refine; i++) {
983ff6701fcSJed Brown       Mat interp_next;
984ff6701fcSJed Brown 
985ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
986ccaff030SJeremy L Thompson       CHKERRQ(ierr);
987ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
988ff6701fcSJed Brown       d = (d + 1) / 2;
989ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
990ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
991ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
992ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
993ff6701fcSJed Brown       if (!i) interpviz = interp_next;
994ff6701fcSJed Brown       else {
995ff6701fcSJed Brown         Mat C;
996ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
997ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
998ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
999ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1000ff6701fcSJed Brown         interpviz = C;
1001ff6701fcSJed Brown       }
1002ff6701fcSJed Brown     }
1003cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1004ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1005cb3e2689Svaleriabarra     }
1006ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1007ccaff030SJeremy L Thompson   }
1008ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1009ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1010ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1011ccaff030SJeremy L Thompson   lnodes /= ncompq;
1012ccaff030SJeremy L Thompson 
10130c6c0b13SLeila Ghaffari   {
10140c6c0b13SLeila Ghaffari     // Print grid information
1015ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1016ccaff030SJeremy L Thompson     int comm_size;
1017ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1018ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1019ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
1020ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1021ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1022ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
10230c6c0b13SLeila Ghaffari     if (!test) {
10240c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d ranks\n",
10250c6c0b13SLeila Ghaffari                          gdofs, odofs, comm_size); CHKERRQ(ierr);
10260c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr);
10270c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str);
1028ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1029ccaff030SJeremy L Thompson     }
1030ccaff030SJeremy L Thompson 
10310c6c0b13SLeila Ghaffari   }
10320c6c0b13SLeila Ghaffari 
1033ccaff030SJeremy L Thompson   // Set up global mass vector
1034ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q,&user->M); CHKERRQ(ierr);
1035ccaff030SJeremy L Thompson 
10360c6c0b13SLeila Ghaffari   // Set up CEED
1037ccaff030SJeremy L Thompson   // CEED Bases
1038ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
1039ea6e0f84SLeila Ghaffari   numP_Vol = degreeVol + 1;
1040ea6e0f84SLeila Ghaffari   numQ_Vol = numP_Vol + qextraVol;
1041ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP_Vol, numQ_Vol, CEED_GAUSS,
1042ea6e0f84SLeila Ghaffari                                   &basisqVol);
1043ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ_Vol, CEED_GAUSS,
1044ea6e0f84SLeila Ghaffari                                   &basisxVol);
1045ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP_Vol,
1046ea6e0f84SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcVol);
1047ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1048ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1049ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1050ccaff030SJeremy L Thompson 
1051ccaff030SJeremy L Thompson   // CEED Restrictions
1052c96c872fSLeila Ghaffari   // Restrictions on the Volume
10536a0edaf9SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dim, 0, 0, 0, numP_Vol, numQ_Vol,
10546a0edaf9SLeila Ghaffari                                  qdatasizeVol, &restrictqVol, &restrictxVol,
10556a0edaf9SLeila Ghaffari                                  &restrictqdiVol); CHKERRQ(ierr);
1056ccaff030SJeremy L Thompson 
1057ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1058ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1059ccaff030SJeremy L Thompson 
1060ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1061bd910870SLeila Ghaffari   CeedInt NqptsVol;
1062c96c872fSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqVol, &NqptsVol);
1063c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqVol, &localNelemVol);
10648b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
1065*cfa64770SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &q0ceed, NULL);
1066ccaff030SJeremy L Thompson 
1067ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1068ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1069ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1070ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1071ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
10728b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1073ccaff030SJeremy L Thompson 
1074ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1075ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1076ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1077ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1078ccaff030SJeremy L Thompson 
1079ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1080ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1081ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1082ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1083ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1084ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
10858b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1086ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1087ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1088ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1089ccaff030SJeremy L Thompson   }
1090ccaff030SJeremy L Thompson 
1091ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1092ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1093ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1094ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1095ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1096ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1097ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
10988b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1099ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1100ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1101ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1102ccaff030SJeremy L Thompson   }
1103ccaff030SJeremy L Thompson 
1104ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1105ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
1106ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictxVol, basisxVol, CEED_VECTOR_ACTIVE);
1107ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
1108ea6e0f84SLeila Ghaffari                        basisxVol, CEED_VECTOR_NONE);
11098b982baeSLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdiVol,
1110ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1111ccaff030SJeremy L Thompson 
1112ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1113ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
1114ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictxVol, basisxcVol, CEED_VECTOR_ACTIVE);
1115ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictqVol,
1116ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1117ccaff030SJeremy L Thompson 
1118ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qceed, NULL);
1119ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qdotceed, NULL);
1120ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->gceed, NULL);
1121ccaff030SJeremy L Thompson 
1122ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1123ccaff030SJeremy L Thompson     CeedOperator op;
1124ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
1125ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1126ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
11278b982baeSLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdiVol,
11288b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
1129ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1130ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1131ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1132ccaff030SJeremy L Thompson     user->op_rhs = op;
1133ccaff030SJeremy L Thompson   }
1134ccaff030SJeremy L Thompson 
1135ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1136ccaff030SJeremy L Thompson     CeedOperator op;
1137ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
1138ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1139ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1140ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictqVol, basisqVol, user->qdotceed);
11418b982baeSLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdiVol,
11428b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
1143ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1144ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1145ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1146ccaff030SJeremy L Thompson     user->op_ifunction = op;
1147ccaff030SJeremy L Thompson   }
1148ccaff030SJeremy L Thompson 
1149*cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1150*cfa64770SLeila Ghaffari   // Outflow Boundary Condition
11516a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
11526a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
1153*cfa64770SLeila Ghaffari   CeedInt numP_Sur, numQ_Sur;
11546a0edaf9SLeila Ghaffari   CeedInt height = 1;
11556a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
11566a0edaf9SLeila Ghaffari   numP_Sur = degreeSur + 1;
11576a0edaf9SLeila Ghaffari   numQ_Sur = numP_Sur + qextraSur;
1158*cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1159*cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1160*cfa64770SLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
1161*cfa64770SLeila Ghaffari   CeedInt NqptsSur;
1162*cfa64770SLeila Ghaffari   PetscInt localNelemSur[6];
1163*cfa64770SLeila Ghaffari   CeedVector qdataSur[6], qdataSur_[6];
1164*cfa64770SLeila Ghaffari   CeedQFunction qf_setupSur, qf_rhsSur, qf_ifunctionSur;
1165*cfa64770SLeila Ghaffari   CeedOperator op_setupSur[6], op_setupSur_[6];
1166*cfa64770SLeila Ghaffari   PetscInt numOutFlow = bc.noutflow;
1167*cfa64770SLeila Ghaffari   DMLabel domainLabel;
1168*cfa64770SLeila Ghaffari 
1169*cfa64770SLeila Ghaffari   // Get Label for the boundaries
1170*cfa64770SLeila Ghaffari   ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
1171*cfa64770SLeila Ghaffari 
1172*cfa64770SLeila Ghaffari   // CEED bases for the boundaries
11736a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
11746a0edaf9SLeila Ghaffari                                   &basisqSur);
11756a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
11766a0edaf9SLeila Ghaffari                                   &basisxSur);
11776a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
11786a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
11796a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
11806a0edaf9SLeila Ghaffari 
1181*cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
11826a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
11836a0edaf9SLeila Ghaffari                               &qf_setupSur);
11846a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
11856a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
11866a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
11876a0edaf9SLeila Ghaffari 
11886a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
1189*cfa64770SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator on the Surface
11906a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
11916a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
11926a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
11938b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "dq", ncompq*dim, CEED_EVAL_GRAD); //assumed volumetric elements
11946a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
11956a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
1196*cfa64770SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdata", qdatasizeVol, CEED_EVAL_NONE); //assumed volumetric elements
11976a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
11986a0edaf9SLeila Ghaffari   }
11996a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
12006a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
12016a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
12026a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
12036a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
12048b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "dq", ncompq*dim, CEED_EVAL_GRAD); //assumed volumetric elements
12056a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
12066a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
1207*cfa64770SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdata", qdatasizeVol, CEED_EVAL_NONE); //assumed volumetric elements
12086a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
12096a0edaf9SLeila Ghaffari   }
1210*cfa64770SLeila Ghaffari   // Create CEED Operator for each face
1211*cfa64770SLeila Ghaffari   for(CeedInt i=0; i<numOutFlow; i++){
1212*cfa64770SLeila Ghaffari 
1213*cfa64770SLeila Ghaffari     ierr = GetRestrictionForDomain(ceed, dm, ncompx, dimSur, height, domainLabel, bc.outflow[i], numP_Sur,
1214*cfa64770SLeila Ghaffari                                    numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
1215*cfa64770SLeila Ghaffari                                    &restrictqdiSur[i]); CHKERRQ(ierr);
1216*cfa64770SLeila Ghaffari     // Create the CEED vectors that will be needed in setup
1217*cfa64770SLeila Ghaffari     CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
1218*cfa64770SLeila Ghaffari     CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
1219*cfa64770SLeila Ghaffari     CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur_[i]);
12206a0edaf9SLeila Ghaffari 
12216a0edaf9SLeila Ghaffari     // Create the operator that builds the quadrature data for the NS operator
1222*cfa64770SLeila Ghaffari     CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
1223*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, CEED_VECTOR_ACTIVE);
1224*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
12256a0edaf9SLeila Ghaffari                          basisxSur, CEED_VECTOR_NONE);
1226*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
12276a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
12286a0edaf9SLeila Ghaffari 
1229*cfa64770SLeila Ghaffari     // Utility operator that builds the quadrature data for computing viscous terms
1230*cfa64770SLeila Ghaffari     CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupSur_[i]);
1231*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "dx", restrictxSur[i], basisxVol, CEED_VECTOR_ACTIVE);
1232*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "weight", CEED_ELEMRESTRICTION_NONE,
1233*cfa64770SLeila Ghaffari                          basisxVol, CEED_VECTOR_NONE);
1234*cfa64770SLeila Ghaffari     CeedOperatorSetField(op_setupSur_[i], "qdata", restrictqdiSur[i],
1235*cfa64770SLeila Ghaffari                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
12366a0edaf9SLeila Ghaffari 
12376a0edaf9SLeila Ghaffari     if (qf_rhsSur) { // Create the RHS physics operator
12386a0edaf9SLeila Ghaffari       CeedOperator op;
12396a0edaf9SLeila Ghaffari       CeedOperatorCreate(ceed, qf_rhsSur, NULL, NULL, &op);
1240*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1241*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "dq", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1242*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdataSur", restrictqdiSur[i],
1243*cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
1244*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "x", restrictxSur[i], basisxSur, xcorners);
1245*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdata", restrictqdiSur[i],
1246*cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur_[i]);
1247*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1248*cfa64770SLeila Ghaffari       user->op_rhs_sur[i] = op;
12496a0edaf9SLeila Ghaffari     }
12506a0edaf9SLeila Ghaffari 
12516a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) { // Create the IFunction operator
12526a0edaf9SLeila Ghaffari       CeedOperator op;
12536a0edaf9SLeila Ghaffari       CeedOperatorCreate(ceed, qf_ifunctionSur, NULL, NULL, &op);
1254*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1255*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "dq", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1256*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdataSur", restrictqdiSur[i],
1257*cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
1258*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "x", restrictxSur[i], basisxSur, xcorners);
1259*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "qdata", restrictqdiSur[i],
1260*cfa64770SLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur_[i]);
1261*cfa64770SLeila Ghaffari       CeedOperatorSetField(op, "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
1262*cfa64770SLeila Ghaffari       user->op_ifunction_sur[i] = op;
1263*cfa64770SLeila Ghaffari     }
12646a0edaf9SLeila Ghaffari   }
12656a0edaf9SLeila Ghaffari   // Composite Operaters
12666a0edaf9SLeila Ghaffari   if (user->op_ifunction_vol) {
1267*cfa64770SLeila Ghaffari     if (numOutFlow>0) {
12686a0edaf9SLeila Ghaffari       // Composite Operators for the IFunction
12696a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_ifunction);
12706a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol);
1271*cfa64770SLeila Ghaffari       for(CeedInt i=0; i<numOutFlow; i++){
1272*cfa64770SLeila Ghaffari         CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_sur[i]);
1273*cfa64770SLeila Ghaffari       }
12746a0edaf9SLeila Ghaffari   } else {
12756a0edaf9SLeila Ghaffari     user->op_ifunction = user->op_ifunction_vol;
12766a0edaf9SLeila Ghaffari     }
12776a0edaf9SLeila Ghaffari   }
12786a0edaf9SLeila Ghaffari   if (user->op_rhs_vol) {
1279*cfa64770SLeila Ghaffari     if (numOutFlow == 1) {
12806a0edaf9SLeila Ghaffari       // Composite Operators for the RHS
12816a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_rhs);
12826a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_vol);
1283*cfa64770SLeila Ghaffari       for(CeedInt i=0; i<numOutFlow; i++){
1284*cfa64770SLeila Ghaffari         CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_sur[i]);
1285*cfa64770SLeila Ghaffari       }
12866a0edaf9SLeila Ghaffari   } else {
12876a0edaf9SLeila Ghaffari     user->op_rhs = user->op_rhs_vol;
12886a0edaf9SLeila Ghaffari     }
12896a0edaf9SLeila Ghaffari   }
1290*cfa64770SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1291ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1292ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1293ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1294ccaff030SJeremy L Thompson     .CtauS = CtauS,
1295ccaff030SJeremy L Thompson     .strong_form = strong_form,
1296ccaff030SJeremy L Thompson     .stabilization = stab,
1297ccaff030SJeremy L Thompson   };
1298ccaff030SJeremy L Thompson   switch (problemChoice) {
1299ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1300ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1301ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1302ccaff030SJeremy L Thompson           sizeof ctxNS);
13036a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
13046a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
13056a0edaf9SLeila Ghaffari           sizeof ctxNS);
1306ccaff030SJeremy L Thompson     break;
1307ccaff030SJeremy L Thompson   case NS_ADVECTION:
1308ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1309ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1310ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1311ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1312ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
13136a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
13146a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
13156a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
13166a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1317ccaff030SJeremy L Thompson   }
1318ccaff030SJeremy L Thompson 
1319ccaff030SJeremy L Thompson   // Set up PETSc context
1320ccaff030SJeremy L Thompson   // Set up units structure
1321ccaff030SJeremy L Thompson   units->meter = meter;
1322ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1323ccaff030SJeremy L Thompson   units->second = second;
1324ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1325ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1326ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1327ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1328ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1329ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1330ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1331ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1332ccaff030SJeremy L Thompson 
1333ccaff030SJeremy L Thompson   // Set up user structure
1334ccaff030SJeremy L Thompson   user->comm = comm;
1335ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1336ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1337ccaff030SJeremy L Thompson   user->units = units;
1338ccaff030SJeremy L Thompson   user->dm = dm;
1339ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1340ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1341ccaff030SJeremy L Thompson   user->ceed = ceed;
1342ccaff030SJeremy L Thompson 
13438b982baeSLeila Ghaffari   // Calculate qdata and ICs
1344ccaff030SJeremy L Thompson   // Set up state global and local vectors
1345ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1346ccaff030SJeremy L Thompson 
1347*cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1348ccaff030SJeremy L Thompson 
1349ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1350ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
13518b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
13528b982baeSLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictqVol, basisqVol, restrictqdiVol, qdata,
1353ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1354ccaff030SJeremy L Thompson 
1355*cfa64770SLeila Ghaffari   ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictqVol,
13560c6c0b13SLeila Ghaffari                                &ctxSetup, 0.0);
1357ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1358ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1359ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1360ccaff030SJeremy L Thompson     // slower execution.
1361ccaff030SJeremy L Thompson     Vec Qbc;
1362ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1363ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1364ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1365ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1366ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1367ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1368ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
13690c6c0b13SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_NS); CHKERRQ(ierr);
1370ccaff030SJeremy L Thompson   }
1371ccaff030SJeremy L Thompson 
1372ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1373ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1374ccaff030SJeremy L Thompson   // Gather initial Q values
1375ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1376ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1377ccaff030SJeremy L Thompson     PetscViewer viewer;
1378ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1379ccaff030SJeremy L Thompson     // Read input
1380ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1381ccaff030SJeremy L Thompson                          user->outputfolder);
1382ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1383ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1384ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1385ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1386ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
13870c6c0b13SLeila Ghaffari   } else {
13880c6c0b13SLeila Ghaffari     //ierr = DMLocalToGlobal(dm, Qloc, INSERT_VALUES, Q);CHKERRQ(ierr);
1389ccaff030SJeremy L Thompson   }
1390ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1391ccaff030SJeremy L Thompson 
1392ccaff030SJeremy L Thompson   // Create and setup TS
1393ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1394ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1395ccaff030SJeremy L Thompson   if (implicit) {
1396ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1397ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1398ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1399ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1400ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1401ccaff030SJeremy L Thompson     }
1402ccaff030SJeremy L Thompson   } else {
1403ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1404ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1405ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1406ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1407ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1408ccaff030SJeremy L Thompson   }
1409ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1410ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1411ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
14120c6c0b13SLeila Ghaffari   if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);}
1413ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1414ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1415ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1416ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1417ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1418ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
14190c6c0b13SLeila Ghaffari     if (!test) {
1420ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1421ccaff030SJeremy L Thompson     }
1422ccaff030SJeremy L Thompson   } else { // continue from time of last output
1423ccaff030SJeremy L Thompson     PetscReal time;
1424ccaff030SJeremy L Thompson     PetscInt count;
1425ccaff030SJeremy L Thompson     PetscViewer viewer;
1426ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1427ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1428ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1429ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1430ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1431ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1432ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1433ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1434ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1435ccaff030SJeremy L Thompson   }
14360c6c0b13SLeila Ghaffari   if (!test) {
1437ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1438ccaff030SJeremy L Thompson   }
1439ccaff030SJeremy L Thompson 
1440ccaff030SJeremy L Thompson   // Solve
1441ccaff030SJeremy L Thompson   start = MPI_Wtime();
1442ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1443ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1444ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1445ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts,&ftime); CHKERRQ(ierr);
1446ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1447ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
14480c6c0b13SLeila Ghaffari   if (!test) {
1449ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
14500c6c0b13SLeila Ghaffari                        "Time taken for solution: %g\n",
1451ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1452ccaff030SJeremy L Thompson   }
1453ccaff030SJeremy L Thompson 
1454ccaff030SJeremy L Thompson   // Get error
14550c6c0b13SLeila Ghaffari   if (problem->non_zero_time && !test) {
1456ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1457ccaff030SJeremy L Thompson     PetscReal norm;
1458ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1459ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1460ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1461ccaff030SJeremy L Thompson 
1462*cfa64770SLeila Ghaffari     ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
1463ea6e0f84SLeila Ghaffari                                  restrictqVol, &ctxSetup, ftime); CHKERRQ(ierr);
1464ccaff030SJeremy L Thompson 
1465ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1466ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1467*cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1468ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1469ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1470ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
1471ccaff030SJeremy L Thompson   }
1472ccaff030SJeremy L Thompson 
1473ccaff030SJeremy L Thompson   // Output Statistics
1474ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
14750c6c0b13SLeila Ghaffari   if (!test) {
1476ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1477ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1478ccaff030SJeremy L Thompson                        steps,(double)ftime); CHKERRQ(ierr);
1479ccaff030SJeremy L Thompson   }
14809cf88b28Svaleriabarra 
1481ccaff030SJeremy L Thompson   // Clean up libCEED
14828b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1483ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1484ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1485ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1486ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
1487ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisqVol);
1488ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxVol);
1489ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxcVol);
1490ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqVol);
1491ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxVol);
1492ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiVol);
1493ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1494ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1495ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1496ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1497ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1498ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
14996a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
15006a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
15016a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
15026a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
15036a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
15046a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
15056a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
15066a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
15076a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
1508*cfa64770SLeila Ghaffari   for(CeedInt i=0; i<numOutFlow; i++){
1509*cfa64770SLeila Ghaffari     CeedVectorDestroy(&qdataSur[i]);
1510*cfa64770SLeila Ghaffari     CeedVectorDestroy(&qdataSur_[i]);
1511*cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictqSur[i]);
1512*cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictxSur[i]);
1513*cfa64770SLeila Ghaffari     CeedElemRestrictionDestroy(&restrictqdiSur[i]);
1514*cfa64770SLeila Ghaffari     CeedOperatorDestroy(&op_setupSur[i]);
1515*cfa64770SLeila Ghaffari     CeedOperatorDestroy(&op_setupSur_[i]);
1516*cfa64770SLeila Ghaffari     CeedOperatorDestroy(&user->op_rhs_sur[i]);
1517*cfa64770SLeila Ghaffari     CeedOperatorDestroy(&user->op_ifunction_sur[i]);
1518*cfa64770SLeila Ghaffari   }
1519ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1520ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1521ccaff030SJeremy L Thompson 
1522ccaff030SJeremy L Thompson   // Clean up PETSc
1523ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1524ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1525ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1526ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1527ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1528ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1529ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1530ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1531ccaff030SJeremy L Thompson   return PetscFinalize();
1532ccaff030SJeremy L Thompson }
1533ccaff030SJeremy L Thompson 
1534