xref: /libCEED/examples/fluids/navierstokes.c (revision 6a0edaf9644a8ecdfc86f4a6d749cd7e1b2dc07e)
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 {
79ea6e0f84SLeila 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,
93c96c872fSLeila Ghaffari     .qdatasizeSur              = 5,
94b0137797SLeila Ghaffari     .setupVol                  = Setup,
95b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
96b0137797SLeila Ghaffari   //.setupSur                  = SetupBoundary,
97b0137797SLeila 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,
102ea6e0f84SLeila Ghaffari   //.applySur_rhs              = DC_Sur,
103ea6e0f84SLeila 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,
114c96c872fSLeila Ghaffari     .qdatasizeSur              = 5,
115b0137797SLeila Ghaffari     .setupVol                  = Setup,
116b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
117b0137797SLeila Ghaffari   //.setupSur                  = SetupBoundary,
118b0137797SLeila 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,
123ea6e0f84SLeila Ghaffari   //.applySur_rhs              = Advection_Sur,
124ea6e0f84SLeila 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,
135b0137797SLeila Ghaffari     .qdatasizeSur              = 1,
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;
168ea6e0f84SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs_sur, op_rhs,
169ea6e0f84SLeila Ghaffari                op_ifunction_vol, op_ifunction_sur, 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_ {
193ccaff030SJeremy L Thompson   PetscInt nwall, nslip[3];
1940c6c0b13SLeila Ghaffari   PetscInt walls[10], slips[3][10];
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,
2040c6c0b13SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt 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);
218ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm,&section); CHKERRQ(ierr);
219ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
220ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
221ccaff030SJeremy L Thompson   fieldoff[0] = 0;
222ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
223ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
224ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
225ccaff030SJeremy L Thompson   }
226ccaff030SJeremy L Thompson 
2270c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2280c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2290c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2300c6c0b13SLeila Ghaffari   if (domainLabel) {
2310c6c0b13SLeila Ghaffari     IS domainIS;
2320c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2330c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2340c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2350c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2360c6c0b13SLeila Ghaffari   } else {
2370c6c0b13SLeila Ghaffari     iterIS = depthIS;
2380c6c0b13SLeila Ghaffari   }
2390c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2400c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
241ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2420c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2430c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
244ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
2450c6c0b13SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
2460c6c0b13SLeila Ghaffari                                    &indices, NULL); CHKERRQ(ierr);
2470c6c0b13SLeila Ghaffari     if (numindices % fieldoff[nfields])
2480c6c0b13SLeila Ghaffari       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
2490c6c0b13SLeila Ghaffari                "Number of closure indices not compatible with Cell %D", c);
250ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
251ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
252ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
253ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
254ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
255ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
256ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
257ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
258ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
259ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
260ccaff030SJeremy L Thompson                      c, i, f, j);
261ccaff030SJeremy L Thompson         }
262ccaff030SJeremy L Thompson       }
263ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
264ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
2656f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
266ccaff030SJeremy L Thompson     }
2670c6c0b13SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
2680c6c0b13SLeila Ghaffari                                        &indices, NULL); CHKERRQ(ierr);
269ccaff030SJeremy L Thompson   }
2700c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
2710c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
2720c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
273ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
2740c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
2750c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
2760c6c0b13SLeila Ghaffari 
277ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
278ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
279ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
2806f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
2816f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
2826f55dfd5Svaleriabarra                             Erestrict);
283ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
284ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
285ccaff030SJeremy L Thompson }
286ccaff030SJeremy L Thompson 
287c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
288bd910870SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt ncompx, CeedInt dim,
289c96c872fSLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value, CeedInt P, CeedInt Q,
290c96c872fSLeila Ghaffari     CeedInt qdatasize, CeedElemRestriction *restrictq,
291c96c872fSLeila Ghaffari     CeedElemRestriction *restrictx, CeedElemRestriction *restrictqdi) {
292c96c872fSLeila Ghaffari 
293c96c872fSLeila Ghaffari   DM dmcoord;
294c96c872fSLeila Ghaffari   CeedInt localNelem;
295c96c872fSLeila Ghaffari   CeedInt Qdim = CeedIntPow(Q, dim);
296c96c872fSLeila Ghaffari   PetscErrorCode ierr;
297c96c872fSLeila Ghaffari 
298c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
299c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
300c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
301c96c872fSLeila Ghaffari   CHKERRQ(ierr);
302c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
303c96c872fSLeila Ghaffari   CHKERRQ(ierr);
304c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
305c96c872fSLeila Ghaffari   CHKERRQ(ierr);
306c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
307c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
308c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
309c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
310c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
311c96c872fSLeila Ghaffari }
312c96c872fSLeila Ghaffari 
313ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
314ccaff030SJeremy L Thompson   PetscErrorCode ierr;
315ccaff030SJeremy L Thompson   PetscInt m;
316ccaff030SJeremy L Thompson 
317ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
318ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
319ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
320ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
321ccaff030SJeremy L Thompson }
322ccaff030SJeremy L Thompson 
323ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
324ccaff030SJeremy L Thompson   PetscErrorCode ierr;
325ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
326ccaff030SJeremy L Thompson   PetscScalar *a;
327ccaff030SJeremy L Thompson 
328ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
329ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
330ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
331ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
3320c6c0b13SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",mpetsc,mceed);
333ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
334ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
335ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
336ccaff030SJeremy L Thompson }
337ccaff030SJeremy L Thompson 
338ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
339ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
340ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
341ccaff030SJeremy L Thompson   PetscErrorCode ierr;
342ccaff030SJeremy L Thompson   Vec Qbc;
343ccaff030SJeremy L Thompson 
344ccaff030SJeremy L Thompson   PetscFunctionBegin;
345ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
346ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
347ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
348ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
349ccaff030SJeremy L Thompson }
350ccaff030SJeremy L Thompson 
351ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
352ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
353ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
354ccaff030SJeremy L Thompson   PetscErrorCode ierr;
355ccaff030SJeremy L Thompson   User user = *(User *)userData;
356ccaff030SJeremy L Thompson   PetscScalar *q, *g;
357ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
358ccaff030SJeremy L Thompson 
359ccaff030SJeremy L Thompson   // Global-to-local
360ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
361ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
362ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
365ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
366ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
367ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
368ccaff030SJeremy L Thompson 
369ccaff030SJeremy L Thompson   // Ceed Vectors
370ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
371ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
373ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
374ccaff030SJeremy L Thompson 
375ccaff030SJeremy L Thompson   // Apply CEED operator
376ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
377ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
378ccaff030SJeremy L Thompson 
379ccaff030SJeremy L Thompson   // Restore vectors
380ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
381ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
382ccaff030SJeremy L Thompson 
383ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
384ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
385ccaff030SJeremy L Thompson 
386ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
387ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
388ccaff030SJeremy L Thompson   CHKERRQ(ierr);
389ccaff030SJeremy L Thompson 
390ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
391ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
392ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
393ccaff030SJeremy L Thompson }
394ccaff030SJeremy L Thompson 
395ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
396ccaff030SJeremy L Thompson                                    void *userData) {
397ccaff030SJeremy L Thompson   PetscErrorCode ierr;
398ccaff030SJeremy L Thompson   User user = *(User *)userData;
399ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
400ccaff030SJeremy L Thompson   PetscScalar *g;
401ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
402ccaff030SJeremy L Thompson 
403ccaff030SJeremy L Thompson   // Global-to-local
404ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
405ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
406ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
408ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
409ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
411ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
412ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
413ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
414ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
415ccaff030SJeremy L Thompson 
416ccaff030SJeremy L Thompson   // Ceed Vectors
417ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
418ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
419ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
420ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
421ccaff030SJeremy L Thompson                      (PetscScalar *)q);
422ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
423ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
424ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
425ccaff030SJeremy L Thompson 
426ccaff030SJeremy L Thompson   // Apply CEED operator
427ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
428ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
429ccaff030SJeremy L Thompson 
430ccaff030SJeremy L Thompson   // Restore vectors
431ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
432ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
433ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson 
435ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
436ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
437ccaff030SJeremy L Thompson 
438ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
439ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
440ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
442ccaff030SJeremy L Thompson }
443ccaff030SJeremy L Thompson 
444ccaff030SJeremy L Thompson // User provided TS Monitor
445ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
446ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
447ccaff030SJeremy L Thompson   User user = ctx;
448ccaff030SJeremy L Thompson   Vec Qloc;
449ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
450ccaff030SJeremy L Thompson   PetscViewer viewer;
451ccaff030SJeremy L Thompson   PetscErrorCode ierr;
452ccaff030SJeremy L Thompson 
453ccaff030SJeremy L Thompson   // Set up output
454ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
455ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
456ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
457ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
458ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
459ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
460ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
461ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
462ccaff030SJeremy L Thompson 
463ccaff030SJeremy L Thompson   // Output
464ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
465ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
466ccaff030SJeremy L Thompson   CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
468ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
469ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
470ccaff030SJeremy L Thompson   if (user->dmviz) {
471ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
472ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
473ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
474ccaff030SJeremy L Thompson 
475ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
476ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
478ccaff030SJeremy L Thompson     CHKERRQ(ierr);
479ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
480ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
482ccaff030SJeremy L Thompson     CHKERRQ(ierr);
483ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
484ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
485ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
486ccaff030SJeremy L Thompson     CHKERRQ(ierr);
487ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
488ccaff030SJeremy L Thompson                               filepath_refined,
489ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
490ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
491ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
494ccaff030SJeremy L Thompson   }
4950c6c0b13SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
496ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson 
498ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
499ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
500ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
501ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
502ccaff030SJeremy L Thompson   CHKERRQ(ierr);
503ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
504ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
505ccaff030SJeremy L Thompson 
506ccaff030SJeremy L Thompson   // Save time stamp
507ccaff030SJeremy L Thompson   // Dimensionalize time back
508ccaff030SJeremy L Thompson   time /= user->units->second;
509ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
510ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
511ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
512ccaff030SJeremy L Thompson   CHKERRQ(ierr);
513ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
514ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
515ccaff030SJeremy L Thompson   #else
516ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
517ccaff030SJeremy L Thompson   #endif
518ccaff030SJeremy L Thompson   CHKERRQ(ierr);
519ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
520ccaff030SJeremy L Thompson 
521ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
522ccaff030SJeremy L Thompson }
523ccaff030SJeremy L Thompson 
5240c6c0b13SLeila Ghaffari static PetscErrorCode ICs_PetscMultiplicity(CeedOperator op_ics,
525ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
526ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
527ccaff030SJeremy L Thompson   PetscErrorCode ierr;
528ccaff030SJeremy L Thompson   CeedVector multlvec;
529ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
530ccaff030SJeremy L Thompson 
531ccaff030SJeremy L Thompson   ctxSetup->time = time;
532ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
535ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
536ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson 
538ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
539ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
541ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
543ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
544ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
545ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
546ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
547ccaff030SJeremy L Thompson   CHKERRQ(ierr);
548ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
549ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
550ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
554ccaff030SJeremy L Thompson }
555ccaff030SJeremy L Thompson 
556ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
557ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
558ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
559ccaff030SJeremy L Thompson   PetscErrorCode ierr;
560ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
561ccaff030SJeremy L Thompson   CeedOperator op_mass;
562ccaff030SJeremy L Thompson   CeedVector mceed;
563ccaff030SJeremy L Thompson   Vec Mloc;
564ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
565ccaff030SJeremy L Thompson 
566ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
567ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
568ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
569ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
570ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
571ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
572ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
573ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
574ccaff030SJeremy L Thompson 
575ccaff030SJeremy L Thompson   // Create the mass operator
576ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
577ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
578ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
579ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
580ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
581ccaff030SJeremy L Thompson 
582ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
584ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
585ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
586ccaff030SJeremy L Thompson 
587ccaff030SJeremy L Thompson   {
588ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
589ccaff030SJeremy L Thompson     CeedVector onesvec;
590ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
591ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
592ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
593ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
594ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
595ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
596ccaff030SJeremy L Thompson   }
597ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
598ccaff030SJeremy L Thompson 
599ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
600ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
601ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
602ccaff030SJeremy L Thompson 
603ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
604ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
605ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
606ccaff030SJeremy L Thompson }
607ccaff030SJeremy L Thompson 
6080c6c0b13SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
609ff6701fcSJed Brown                        SimpleBC bc, void *ctxSetup) {
610ccaff030SJeremy L Thompson   PetscErrorCode ierr;
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
613ccaff030SJeremy L Thompson   {
614ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
615ccaff030SJeremy L Thompson     PetscFE fe;
6160c6c0b13SLeila Ghaffari     PetscSpace fespace;
617ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
618ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
619ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
62032ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
621ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
622ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
623ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
6240c6c0b13SLeila Ghaffari     /* Wall boundary conditions are zero velocity and zero flux for density and energy */
6250c6c0b13SLeila Ghaffari     {
6260c6c0b13SLeila Ghaffari       PetscInt comps[3] = {1, 2, 3};
6270c6c0b13SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
6280c6c0b13SLeila Ghaffari                            3, comps, (void(*)(void))problem->bc,
6290c6c0b13SLeila Ghaffari                            bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
6300c6c0b13SLeila Ghaffari     }
63107af6069Svaleriabarra     {
63207af6069Svaleriabarra       PetscInt comps[1] = {1};
63307af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
63407af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
63507af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
63607af6069Svaleriabarra       comps[0] = 2;
63707af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
63807af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
63907af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
64007af6069Svaleriabarra       comps[0] = 3;
64107af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
64207af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
64307af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
64407af6069Svaleriabarra     }
645ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE,NULL);
646ccaff030SJeremy L Thompson     CHKERRQ(ierr);
6470c6c0b13SLeila Ghaffari     ierr = PetscFEGetBasisSpace(fe, &fespace); CHKERRQ(ierr);
648ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
649ccaff030SJeremy L Thompson   }
650ccaff030SJeremy L Thompson   {
651ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
652ccaff030SJeremy L Thompson     PetscSection section;
653ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
654ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
656ccaff030SJeremy L Thompson     CHKERRQ(ierr);
657ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
658ccaff030SJeremy L Thompson     CHKERRQ(ierr);
659ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
660ccaff030SJeremy L Thompson     CHKERRQ(ierr);
661ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
662ccaff030SJeremy L Thompson     CHKERRQ(ierr);
663ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
664ccaff030SJeremy L Thompson     CHKERRQ(ierr);
665ccaff030SJeremy L Thompson   }
666ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
667ccaff030SJeremy L Thompson }
668ccaff030SJeremy L Thompson 
669ccaff030SJeremy L Thompson int main(int argc, char **argv) {
670ccaff030SJeremy L Thompson   PetscInt ierr;
671ccaff030SJeremy L Thompson   MPI_Comm comm;
6720c6c0b13SLeila Ghaffari   DM dm, dmcoord, dmviz, dmvizfine;
673ccaff030SJeremy L Thompson   Mat interpviz;
674ccaff030SJeremy L Thompson   TS ts;
675ccaff030SJeremy L Thompson   TSAdapt adapt;
676ccaff030SJeremy L Thompson   User user;
677ccaff030SJeremy L Thompson   Units units;
678ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
679c96c872fSLeila Ghaffari   PetscInt localNelemVol, localNelemSur, lnodes, steps;
680ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
681ccaff030SJeremy L Thompson   PetscMPIInt rank;
682ccaff030SJeremy L Thompson   PetscScalar ftime;
683ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
684ccaff030SJeremy L Thompson   Ceed ceed;
685ea6e0f84SLeila Ghaffari   CeedInt numP_Vol, numP_Sur, numQ_Vol, numQ_Sur;
686c96c872fSLeila Ghaffari   CeedVector xcorners, qdataVol, qdataSur, q0ceedVol, q0ceedSur;
687ea6e0f84SLeila Ghaffari   CeedBasis basisxVol, basisxcVol, basisqVol, basisxSur, basisxcSur, basisqSur;
688bd910870SLeila Ghaffari   CeedElemRestriction restrictxVol, restrictqVol, restrictqdiVol,
689bd910870SLeila Ghaffari                       restrictxSur, restrictqSur, restrictqdiSur;
690ea6e0f84SLeila Ghaffari   CeedQFunction qf_setupVol, qf_setupSur, qf_ics, qf_rhsVol, qf_rhsSur,
691ea6e0f84SLeila Ghaffari                 qf_ifunctionVol, qf_ifunctionSur;
692ea6e0f84SLeila Ghaffari   CeedOperator op_setupVol, op_setupSur, op_ics;
693ccaff030SJeremy L Thompson   CeedScalar Rd;
694ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
695ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
696ccaff030SJeremy L Thompson   problemType problemChoice;
697ccaff030SJeremy L Thompson   problemData *problem = NULL;
698ccaff030SJeremy L Thompson   StabilizationType stab;
6990c6c0b13SLeila Ghaffari   PetscBool   test, implicit;
700cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
701ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
7020c6c0b13SLeila Ghaffari     .nwall = 6,
7030c6c0b13SLeila Ghaffari     .walls = {1,2,3,4,5,6},
704ccaff030SJeremy L Thompson   };
705ccaff030SJeremy L Thompson   double start, cpu_time_used;
706ccaff030SJeremy L Thompson 
707ccaff030SJeremy L Thompson   // Create the libCEED contexts
708ccaff030SJeremy L Thompson   PetscScalar meter     = 1e-2;     // 1 meter in scaled length units
709ccaff030SJeremy L Thompson   PetscScalar second    = 1e-2;     // 1 second in scaled time units
710ccaff030SJeremy L Thompson   PetscScalar kilogram  = 1e-6;     // 1 kilogram in scaled mass units
711ccaff030SJeremy L Thompson   PetscScalar Kelvin    = 1;        // 1 Kelvin in scaled temperature units
712ccaff030SJeremy L Thompson   CeedScalar theta0     = 300.;     // K
713ccaff030SJeremy L Thompson   CeedScalar thetaC     = -15.;     // K
714ccaff030SJeremy L Thompson   CeedScalar P0         = 1.e5;     // Pa
715ccaff030SJeremy L Thompson   CeedScalar N          = 0.01;     // 1/s
716ccaff030SJeremy L Thompson   CeedScalar cv         = 717.;     // J/(kg K)
717ccaff030SJeremy L Thompson   CeedScalar cp         = 1004.;    // J/(kg K)
718ccaff030SJeremy L Thompson   CeedScalar g          = 9.81;     // m/s^2
719ccaff030SJeremy L Thompson   CeedScalar lambda     = -2./3.;   // -
720ccaff030SJeremy L Thompson   CeedScalar mu         = 75.;      // Pa s, dynamic viscosity
721ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
722ccaff030SJeremy L Thompson   CeedScalar k          = 0.02638;  // W/(m K)
723ccaff030SJeremy L Thompson   CeedScalar CtauS      = 0.;       // dimensionless
724ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
725ccaff030SJeremy L Thompson   PetscScalar lx        = 8000.;    // m
726ccaff030SJeremy L Thompson   PetscScalar ly        = 8000.;    // m
727ccaff030SJeremy L Thompson   PetscScalar lz        = 4000.;    // m
728ccaff030SJeremy L Thompson   CeedScalar rc         = 1000.;    // m (Radius of bubble)
729ccaff030SJeremy L Thompson   PetscScalar resx      = 1000.;    // m (resolution in x)
730ccaff030SJeremy L Thompson   PetscScalar resy      = 1000.;    // m (resolution in y)
731ccaff030SJeremy L Thompson   PetscScalar resz      = 1000.;    // m (resolution in z)
732ccaff030SJeremy L Thompson   PetscInt outputfreq   = 10;       // -
733ccaff030SJeremy L Thompson   PetscInt contsteps    = 0;        // -
734ea6e0f84SLeila Ghaffari   PetscInt degreeVol    = 1;        // -
735ea6e0f84SLeila Ghaffari   PetscInt degreeSur    = 1;        // -
736ea6e0f84SLeila Ghaffari   PetscInt qextraVol    = 2;        // -
737ea6e0f84SLeila Ghaffari   PetscInt qextraSur    = 2;        // -
7380c6c0b13SLeila Ghaffari   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
7390c6c0b13SLeila Ghaffari                                   DM_BOUNDARY_NONE
7400c6c0b13SLeila Ghaffari                                  };
741ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
742ccaff030SJeremy L Thompson 
743ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
744ccaff030SJeremy L Thompson   if (ierr) return ierr;
745ccaff030SJeremy L Thompson 
746ccaff030SJeremy L Thompson   // Allocate PETSc context
747ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
748ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
749ccaff030SJeremy L Thompson 
750ccaff030SJeremy L Thompson   // Parse command line options
751ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
752ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
753ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
754ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
755ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
756ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
7570c6c0b13SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
7580c6c0b13SLeila Ghaffari                           NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr);
759ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
760ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
761ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
762ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
763ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
764ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
765ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
766ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
767ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
768ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
769ccaff030SJeremy L Thompson   CHKERRQ(ierr);
770ccaff030SJeremy L Thompson   {
771ccaff030SJeremy L Thompson     PetscInt len;
772ccaff030SJeremy L Thompson     PetscBool flg;
773ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
774ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
775ccaff030SJeremy L Thompson                                 NULL, bc.walls,
776ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
777ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
778ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
779ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
780ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
781ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
782ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
783ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
784ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
785ccaff030SJeremy L Thompson                                    &len), &flg);
786ccaff030SJeremy L Thompson       CHKERRQ(ierr);
7870c6c0b13SLeila Ghaffari       if (flg) bc.nslip[j] = len;
788ccaff030SJeremy L Thompson     }
789ccaff030SJeremy L Thompson   }
790cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
791cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
792cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
793ccaff030SJeremy L Thompson   CHKERRQ(ierr);
794ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
795ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
796ccaff030SJeremy L Thompson   meter = fabs(meter);
797ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
798ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
799ccaff030SJeremy L Thompson   second = fabs(second);
800ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
801ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
802ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
803ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
804ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
805ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
806ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
807ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
808ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
809ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
810ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
811ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
812ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
813ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
814ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
815ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
816ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
817ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
818ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
819ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
820ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
821ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
822ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
823ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
824ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
825ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
826ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
827ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
828ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
829ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
830ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
831ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
832ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
833ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
834ccaff030SJeremy L Thompson   CHKERRQ(ierr);
835ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
836ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
837ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
838ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
839ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
840ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
841ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
842ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
843ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
844ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
845ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
846ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
847ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
848ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
849ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
8500c6c0b13SLeila Ghaffari   ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction",
8510c6c0b13SLeila Ghaffari                                NULL, DMBoundaryTypes, (PetscEnum *)periodicity,
8520c6c0b13SLeila Ghaffari                                &n, NULL); CHKERRQ(ierr);
8530c6c0b13SLeila Ghaffari   n = problem->dim;
854ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
855ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
856ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
857ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
858ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
859ccaff030SJeremy L Thompson   n = problem->dim;
860ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
861ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
862ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
863ccaff030SJeremy L Thompson   {
864ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
865ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
866ccaff030SJeremy L Thompson     if (norm > 0) {
867ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
868ccaff030SJeremy L Thompson     }
869ccaff030SJeremy L Thompson   }
870ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
871ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
872ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
873ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
874ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
875ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-degree_volume", "Polynomial degree of finite elements for the Volume",
876ea6e0f84SLeila Ghaffari                          NULL, degreeVol, &degreeVol, NULL); CHKERRQ(ierr);
877ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-qextra_volume", "Number of extra quadrature points for the Volume",
878ea6e0f84SLeila Ghaffari                          NULL, qextraVol, &qextraVol, NULL); CHKERRQ(ierr);
8790c6c0b13SLeila Ghaffari   PetscStrncpy(user->outputfolder, ".", 2);
880ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
881ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
882ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
883ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
884ccaff030SJeremy L Thompson 
885ccaff030SJeremy L Thompson   // Define derived units
886ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
887ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
888ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
889ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
890ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
891ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
892ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
893ccaff030SJeremy L Thompson 
894ccaff030SJeremy L Thompson   // Scale variables to desired units
895ccaff030SJeremy L Thompson   theta0 *= Kelvin;
896ccaff030SJeremy L Thompson   thetaC *= Kelvin;
897ccaff030SJeremy L Thompson   P0 *= Pascal;
898ccaff030SJeremy L Thompson   N *= (1./second);
899ccaff030SJeremy L Thompson   cv *= JperkgK;
900ccaff030SJeremy L Thompson   cp *= JperkgK;
901ccaff030SJeremy L Thompson   Rd = cp - cv;
902ccaff030SJeremy L Thompson   g *= mpersquareds;
903ccaff030SJeremy L Thompson   mu *= Pascal * second;
904ccaff030SJeremy L Thompson   k *= WpermK;
905ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
906ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
907ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
908ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
909ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
910ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
911ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
912ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
913ccaff030SJeremy L Thompson 
914ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
915c96c872fSLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol,
916c96c872fSLeila Ghaffari                 qdatasizeSur = problem->qdatasizeSur;
917ccaff030SJeremy L Thompson   // Set up the libCEED context
918ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup =  {
919ccaff030SJeremy L Thompson     .theta0 = theta0,
920ccaff030SJeremy L Thompson     .thetaC = thetaC,
921ccaff030SJeremy L Thompson     .P0 = P0,
922ccaff030SJeremy L Thompson     .N = N,
923ccaff030SJeremy L Thompson     .cv = cv,
924ccaff030SJeremy L Thompson     .cp = cp,
925ccaff030SJeremy L Thompson     .Rd = Rd,
926ccaff030SJeremy L Thompson     .g = g,
927ccaff030SJeremy L Thompson     .rc = rc,
928ccaff030SJeremy L Thompson     .lx = lx,
929ccaff030SJeremy L Thompson     .ly = ly,
930ccaff030SJeremy L Thompson     .lz = lz,
9310c6c0b13SLeila Ghaffari     .periodicity0 = periodicity[0],
9320c6c0b13SLeila Ghaffari     .periodicity1 = periodicity[1],
9330c6c0b13SLeila Ghaffari     .periodicity2 = periodicity[2],
934ccaff030SJeremy L Thompson     .center[0] = center[0],
935ccaff030SJeremy L Thompson     .center[1] = center[1],
936ccaff030SJeremy L Thompson     .center[2] = center[2],
937ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
938ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
939ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
940ccaff030SJeremy L Thompson     .time = 0,
941ccaff030SJeremy L Thompson   };
942ccaff030SJeremy L Thompson 
943ccaff030SJeremy L Thompson   {
944ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
945ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
9460c6c0b13SLeila Ghaffari                                periodicity, PETSC_TRUE, &dm);
947ccaff030SJeremy L Thompson     CHKERRQ(ierr);
948ccaff030SJeremy L Thompson   }
9490c6c0b13SLeila Ghaffari   if (1) {
950ccaff030SJeremy L Thompson     DM               dmDist = NULL;
951ccaff030SJeremy L Thompson     PetscPartitioner part;
952ccaff030SJeremy L Thompson 
953ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
954ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
955ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
956ccaff030SJeremy L Thompson     if (dmDist) {
957ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
958ccaff030SJeremy L Thompson       dm  = dmDist;
959ccaff030SJeremy L Thompson     }
960ccaff030SJeremy L Thompson   }
961ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
962ccaff030SJeremy L Thompson 
963ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
964ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
965ea6e0f84SLeila Ghaffari   ierr = SetUpDM(dm, problem, degreeVol, &bc, &ctxSetup); CHKERRQ(ierr);
9660c6c0b13SLeila Ghaffari   if (!test) {
9670c6c0b13SLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
968ea6e0f84SLeila Ghaffari                        "Degree of Volumetric FEM Space: %D\n",
969ea6e0f84SLeila Ghaffari                        degreeVol); CHKERRQ(ierr);
9700c6c0b13SLeila Ghaffari   }
971ccaff030SJeremy L Thompson   dmviz = NULL;
972ccaff030SJeremy L Thompson   interpviz = NULL;
973ccaff030SJeremy L Thompson   if (viz_refine) {
974ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
975ff6701fcSJed Brown 
976ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
977ff6701fcSJed Brown     dmhierarchy[0] = dm;
978ea6e0f84SLeila Ghaffari     for (PetscInt i = 0, d = degreeVol; i < viz_refine; i++) {
979ff6701fcSJed Brown       Mat interp_next;
980ff6701fcSJed Brown 
981ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
982ccaff030SJeremy L Thompson       CHKERRQ(ierr);
983ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
984ff6701fcSJed Brown       d = (d + 1) / 2;
985ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
986ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
987ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
988ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
989ff6701fcSJed Brown       if (!i) interpviz = interp_next;
990ff6701fcSJed Brown       else {
991ff6701fcSJed Brown         Mat C;
992ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
993ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
994ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
995ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
996ff6701fcSJed Brown         interpviz = C;
997ff6701fcSJed Brown       }
998ff6701fcSJed Brown     }
999cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1000ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1001cb3e2689Svaleriabarra     }
1002ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1003ccaff030SJeremy L Thompson   }
1004ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1005ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1006ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1007ccaff030SJeremy L Thompson   lnodes /= ncompq;
1008ccaff030SJeremy L Thompson 
10090c6c0b13SLeila Ghaffari   {
10100c6c0b13SLeila Ghaffari     // Print grid information
1011ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1012ccaff030SJeremy L Thompson     int comm_size;
1013ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1014ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1015ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
1016ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1017ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1018ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
10190c6c0b13SLeila Ghaffari     if (!test) {
10200c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d ranks\n",
10210c6c0b13SLeila Ghaffari                          gdofs, odofs, comm_size); CHKERRQ(ierr);
10220c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr);
10230c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str);
1024ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1025ccaff030SJeremy L Thompson     }
1026ccaff030SJeremy L Thompson 
10270c6c0b13SLeila Ghaffari   }
10280c6c0b13SLeila Ghaffari 
1029ccaff030SJeremy L Thompson   // Set up global mass vector
1030ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q,&user->M); CHKERRQ(ierr);
1031ccaff030SJeremy L Thompson 
10320c6c0b13SLeila Ghaffari   // Set up CEED
1033ccaff030SJeremy L Thompson   // CEED Bases
1034ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
1035ea6e0f84SLeila Ghaffari   numP_Vol = degreeVol + 1;
1036ea6e0f84SLeila Ghaffari   numQ_Vol = numP_Vol + qextraVol;
1037ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP_Vol, numQ_Vol, CEED_GAUSS,
1038ea6e0f84SLeila Ghaffari                                   &basisqVol);
1039ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ_Vol, CEED_GAUSS,
1040ea6e0f84SLeila Ghaffari                                   &basisxVol);
1041ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP_Vol,
1042ea6e0f84SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcVol);
1043ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1044ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1045ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1046ccaff030SJeremy L Thompson 
1047ccaff030SJeremy L Thompson   // CEED Restrictions
1048c96c872fSLeila Ghaffari   // Restrictions on the Volume
1049*6a0edaf9SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dim, 0, 0, 0, numP_Vol, numQ_Vol,
1050*6a0edaf9SLeila Ghaffari                                  qdatasizeVol, &restrictqVol, &restrictxVol,
1051*6a0edaf9SLeila Ghaffari                                  &restrictqdiVol); CHKERRQ(ierr);
1052ccaff030SJeremy L Thompson 
1053ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1054ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1055ccaff030SJeremy L Thompson 
1056ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1057bd910870SLeila Ghaffari   CeedInt NqptsVol;
1058c96c872fSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqVol, &NqptsVol);
1059c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqVol, &localNelemVol);
1060c96c872fSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdataVol);
1061c96c872fSLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &q0ceedVol, NULL);
1062ccaff030SJeremy L Thompson 
1063ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1064ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1065ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1066ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1067ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
1068ea6e0f84SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1069ccaff030SJeremy L Thompson 
1070ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1071ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1072ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1073ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1074ccaff030SJeremy L Thompson 
1075ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1076ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1077ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1078ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1079ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1080ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1081ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1082ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1083ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1084ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1085ccaff030SJeremy L Thompson   }
1086ccaff030SJeremy L Thompson 
1087ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1088ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1089ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1090ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1091ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1092ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1093ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
1094ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1095ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1096ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1097ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1098ccaff030SJeremy L Thompson   }
1099ccaff030SJeremy L Thompson 
1100ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1101ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
1102ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictxVol, basisxVol, CEED_VECTOR_ACTIVE);
1103ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
1104ea6e0f84SLeila Ghaffari                        basisxVol, CEED_VECTOR_NONE);
1105ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdataVol", restrictqdiVol,
1106ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1107ccaff030SJeremy L Thompson 
1108ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1109ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
1110ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictxVol, basisxcVol, CEED_VECTOR_ACTIVE);
1111ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictqVol,
1112ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1113ccaff030SJeremy L Thompson 
1114ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qceed, NULL);
1115ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qdotceed, NULL);
1116ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->gceed, NULL);
1117ccaff030SJeremy L Thompson 
1118ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1119ccaff030SJeremy L Thompson     CeedOperator op;
1120ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
1121ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1122ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1123ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1124ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1125ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1126ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1127ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1128ccaff030SJeremy L Thompson     user->op_rhs = op;
1129ccaff030SJeremy L Thompson   }
1130ccaff030SJeremy L Thompson 
1131ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1132ccaff030SJeremy L Thompson     CeedOperator op;
1133ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
1134ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1135ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1136ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictqVol, basisqVol, user->qdotceed);
1137ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1138ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1139ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1140ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1141ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1142ccaff030SJeremy L Thompson     user->op_ifunction = op;
1143ccaff030SJeremy L Thompson   }
1144ccaff030SJeremy L Thompson 
1145*6a0edaf9SLeila Ghaffari   //**************************************************************************************//
1146*6a0edaf9SLeila Ghaffari   // Add boundary Integral (TODO: create a function for all faces)
1147*6a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
1148*6a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
1149*6a0edaf9SLeila Ghaffari   // CEED bases
1150*6a0edaf9SLeila Ghaffari   CeedInt height = 1;
1151*6a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
1152*6a0edaf9SLeila Ghaffari   numP_Sur = degreeSur + 1;
1153*6a0edaf9SLeila Ghaffari   numQ_Sur = numP_Sur + qextraSur;
1154*6a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
1155*6a0edaf9SLeila Ghaffari                                   &basisqSur);
1156*6a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
1157*6a0edaf9SLeila Ghaffari                                   &basisxSur);
1158*6a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
1159*6a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
1160*6a0edaf9SLeila Ghaffari   // CEED Restrictions
1161*6a0edaf9SLeila Ghaffari   // Restriction on one face
1162*6a0edaf9SLeila Ghaffari   DMLabel domainLabel;
1163*6a0edaf9SLeila Ghaffari   ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
1164*6a0edaf9SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dimSur, height, domainLabel, 2, numP_Sur,
1165*6a0edaf9SLeila Ghaffari                                  numQ_Sur, qdatasizeSur, &restrictqSur, &restrictxSur,
1166*6a0edaf9SLeila Ghaffari                                  &restrictqdiSur); CHKERRQ(ierr);
1167*6a0edaf9SLeila Ghaffari 
1168*6a0edaf9SLeila Ghaffari   // Create the CEED vectors that will be needed in setup
1169*6a0edaf9SLeila Ghaffari   CeedInt NqptsSur;
1170*6a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
1171*6a0edaf9SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqSur, &localNelemSur);
1172*6a0edaf9SLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeSur*localNelemSur*NqptsSur, &qdataSur);
1173*6a0edaf9SLeila Ghaffari 
1174*6a0edaf9SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the NS operator
1175*6a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
1176*6a0edaf9SLeila Ghaffari                               &qf_setupSur);
1177*6a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
1178*6a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
1179*6a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
1180*6a0edaf9SLeila Ghaffari 
1181*6a0edaf9SLeila Ghaffari 
1182*6a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
1183*6a0edaf9SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator
1184*6a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
1185*6a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
1186*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
1187*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "dq", ncompq*dimSur, CEED_EVAL_GRAD);
1188*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
1189*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
1190*6a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
1191*6a0edaf9SLeila Ghaffari   }
1192*6a0edaf9SLeila Ghaffari 
1193*6a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
1194*6a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
1195*6a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
1196*6a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
1197*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
1198*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "dq", ncompq*dimSur, CEED_EVAL_GRAD);
1199*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
1200*6a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
1201*6a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
1202*6a0edaf9SLeila Ghaffari   }
1203*6a0edaf9SLeila Ghaffari 
1204*6a0edaf9SLeila Ghaffari   // Create the operator that builds the quadrature data for the NS operator
1205*6a0edaf9SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur);
1206*6a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "dx", restrictxSur, basisxSur, CEED_VECTOR_ACTIVE);
1207*6a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "weight", CEED_ELEMRESTRICTION_NONE,
1208*6a0edaf9SLeila Ghaffari                        basisxSur, CEED_VECTOR_NONE);
1209*6a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "qdataSur", restrictqdiSur,
1210*6a0edaf9SLeila Ghaffari                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1211*6a0edaf9SLeila Ghaffari 
1212*6a0edaf9SLeila Ghaffari 
1213*6a0edaf9SLeila Ghaffari   if (qf_rhsSur) { // Create the RHS physics operator
1214*6a0edaf9SLeila Ghaffari     CeedOperator op;
1215*6a0edaf9SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsSur, NULL, NULL, &op);
1216*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1217*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1218*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "qdataSur", restrictqdiSur,
1219*6a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataSur);
1220*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxSur, basisxSur, xcorners);
1221*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1222*6a0edaf9SLeila Ghaffari     user->op_rhs_sur = op;
1223*6a0edaf9SLeila Ghaffari   }
1224*6a0edaf9SLeila Ghaffari 
1225*6a0edaf9SLeila Ghaffari   if (qf_ifunctionSur) { // Create the IFunction operator
1226*6a0edaf9SLeila Ghaffari     CeedOperator op;
1227*6a0edaf9SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionSur, NULL, NULL, &op);
1228*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1229*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1230*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "qdataSur", restrictqdiSur,
1231*6a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataSur);
1232*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxSur, basisxSur, xcorners);
1233*6a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
1234*6a0edaf9SLeila Ghaffari     user->op_ifunction_sur = op;
1235*6a0edaf9SLeila Ghaffari   }
1236*6a0edaf9SLeila Ghaffari   // Composite Operaters
1237*6a0edaf9SLeila Ghaffari   if (user->op_ifunction_vol) {
1238*6a0edaf9SLeila Ghaffari     if (user->op_ifunction_sur) {
1239*6a0edaf9SLeila Ghaffari       // Composite Operators for the IFunction
1240*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_ifunction);
1241*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol);
1242*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_sur);
1243*6a0edaf9SLeila Ghaffari   } else {
1244*6a0edaf9SLeila Ghaffari     user->op_ifunction = user->op_ifunction_vol;
1245*6a0edaf9SLeila Ghaffari     }
1246*6a0edaf9SLeila Ghaffari   }
1247*6a0edaf9SLeila Ghaffari   if (user->op_rhs_vol) {
1248*6a0edaf9SLeila Ghaffari     if (user->op_rhs_sur) {
1249*6a0edaf9SLeila Ghaffari       // Composite Operators for the RHS
1250*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_rhs);
1251*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_vol);
1252*6a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_sur);
1253*6a0edaf9SLeila Ghaffari   } else {
1254*6a0edaf9SLeila Ghaffari     user->op_rhs = user->op_rhs_vol;
1255*6a0edaf9SLeila Ghaffari     }
1256*6a0edaf9SLeila Ghaffari   }
1257*6a0edaf9SLeila Ghaffari   //**************************************************************************************//
1258ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1259ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1260ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1261ccaff030SJeremy L Thompson     .CtauS = CtauS,
1262ccaff030SJeremy L Thompson     .strong_form = strong_form,
1263ccaff030SJeremy L Thompson     .stabilization = stab,
1264ccaff030SJeremy L Thompson   };
1265ccaff030SJeremy L Thompson   switch (problemChoice) {
1266ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1267ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1268ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1269ccaff030SJeremy L Thompson           sizeof ctxNS);
1270*6a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
1271*6a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
1272*6a0edaf9SLeila Ghaffari           sizeof ctxNS);
1273ccaff030SJeremy L Thompson     break;
1274ccaff030SJeremy L Thompson   case NS_ADVECTION:
1275ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1276ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1277ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1278ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1279ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1280*6a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
1281*6a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
1282*6a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
1283*6a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1284ccaff030SJeremy L Thompson   }
1285ccaff030SJeremy L Thompson 
1286ccaff030SJeremy L Thompson   // Set up PETSc context
1287ccaff030SJeremy L Thompson   // Set up units structure
1288ccaff030SJeremy L Thompson   units->meter = meter;
1289ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1290ccaff030SJeremy L Thompson   units->second = second;
1291ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1292ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1293ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1294ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1295ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1296ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1297ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1298ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1299ccaff030SJeremy L Thompson 
1300ccaff030SJeremy L Thompson   // Set up user structure
1301ccaff030SJeremy L Thompson   user->comm = comm;
1302ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1303ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1304ccaff030SJeremy L Thompson   user->units = units;
1305ccaff030SJeremy L Thompson   user->dm = dm;
1306ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1307ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1308ccaff030SJeremy L Thompson   user->ceed = ceed;
1309ccaff030SJeremy L Thompson 
1310ea6e0f84SLeila Ghaffari   // Calculate qdataVol and ICs
1311ccaff030SJeremy L Thompson   // Set up state global and local vectors
1312ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1313ccaff030SJeremy L Thompson 
1314c96c872fSLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceedVol, Qloc); CHKERRQ(ierr);
1315ccaff030SJeremy L Thompson 
1316ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1317ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
1318ea6e0f84SLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdataVol, CEED_REQUEST_IMMEDIATE);
1319ea6e0f84SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictqVol, basisqVol, restrictqdiVol, qdataVol,
1320ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1321ccaff030SJeremy L Thompson 
1322c96c872fSLeila Ghaffari   ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qloc, Q, restrictqVol,
13230c6c0b13SLeila Ghaffari                                &ctxSetup, 0.0);
1324ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1325ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1326ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1327ccaff030SJeremy L Thompson     // slower execution.
1328ccaff030SJeremy L Thompson     Vec Qbc;
1329ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1330ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1331ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1332ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1333ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1334ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1335ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
13360c6c0b13SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_NS); CHKERRQ(ierr);
1337ccaff030SJeremy L Thompson   }
1338ccaff030SJeremy L Thompson 
1339ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1340ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1341ccaff030SJeremy L Thompson   // Gather initial Q values
1342ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1343ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1344ccaff030SJeremy L Thompson     PetscViewer viewer;
1345ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1346ccaff030SJeremy L Thompson     // Read input
1347ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1348ccaff030SJeremy L Thompson                          user->outputfolder);
1349ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1350ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1351ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1352ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1353ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
13540c6c0b13SLeila Ghaffari   } else {
13550c6c0b13SLeila Ghaffari     //ierr = DMLocalToGlobal(dm, Qloc, INSERT_VALUES, Q);CHKERRQ(ierr);
1356ccaff030SJeremy L Thompson   }
1357ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1358ccaff030SJeremy L Thompson 
1359ccaff030SJeremy L Thompson   // Create and setup TS
1360ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1361ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1362ccaff030SJeremy L Thompson   if (implicit) {
1363ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1364ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1365ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1366ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1367ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1368ccaff030SJeremy L Thompson     }
1369ccaff030SJeremy L Thompson   } else {
1370ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1371ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1372ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1373ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1374ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1375ccaff030SJeremy L Thompson   }
1376ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1377ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1378ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
13790c6c0b13SLeila Ghaffari   if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);}
1380ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1381ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1382ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1383ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1384ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1385ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
13860c6c0b13SLeila Ghaffari     if (!test) {
1387ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1388ccaff030SJeremy L Thompson     }
1389ccaff030SJeremy L Thompson   } else { // continue from time of last output
1390ccaff030SJeremy L Thompson     PetscReal time;
1391ccaff030SJeremy L Thompson     PetscInt count;
1392ccaff030SJeremy L Thompson     PetscViewer viewer;
1393ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1394ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1395ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1396ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1397ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1398ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1399ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1400ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1401ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1402ccaff030SJeremy L Thompson   }
14030c6c0b13SLeila Ghaffari   if (!test) {
1404ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1405ccaff030SJeremy L Thompson   }
1406ccaff030SJeremy L Thompson 
1407ccaff030SJeremy L Thompson   // Solve
1408ccaff030SJeremy L Thompson   start = MPI_Wtime();
1409ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1410ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1411ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1412ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts,&ftime); CHKERRQ(ierr);
1413ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1414ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
14150c6c0b13SLeila Ghaffari   if (!test) {
1416ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
14170c6c0b13SLeila Ghaffari                        "Time taken for solution: %g\n",
1418ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1419ccaff030SJeremy L Thompson   }
1420ccaff030SJeremy L Thompson 
1421ccaff030SJeremy L Thompson   // Get error
14220c6c0b13SLeila Ghaffari   if (problem->non_zero_time && !test) {
1423ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1424ccaff030SJeremy L Thompson     PetscReal norm;
1425ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1426ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1427ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1428ccaff030SJeremy L Thompson 
1429c96c872fSLeila Ghaffari     ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qexactloc, Qexact,
1430ea6e0f84SLeila Ghaffari                                  restrictqVol, &ctxSetup, ftime); CHKERRQ(ierr);
1431ccaff030SJeremy L Thompson 
1432ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1433ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1434c96c872fSLeila Ghaffari     CeedVectorDestroy(&q0ceedVol);
1435ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1436ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1437ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
1438ccaff030SJeremy L Thompson   }
1439ccaff030SJeremy L Thompson 
1440ccaff030SJeremy L Thompson   // Output Statistics
1441ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
14420c6c0b13SLeila Ghaffari   if (!test) {
1443ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1444ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1445ccaff030SJeremy L Thompson                        steps,(double)ftime); CHKERRQ(ierr);
1446ccaff030SJeremy L Thompson   }
14479cf88b28Svaleriabarra 
1448ccaff030SJeremy L Thompson   // Clean up libCEED
1449ea6e0f84SLeila Ghaffari   CeedVectorDestroy(&qdataVol);
1450ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1451ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1452ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1453ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
1454ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisqVol);
1455ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxVol);
1456ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxcVol);
1457ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqVol);
1458ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxVol);
1459ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiVol);
1460ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1461ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1462ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1463ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1464ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1465ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
1466*6a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
1467*6a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
1468*6a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
1469*6a0edaf9SLeila Ghaffari 
1470*6a0edaf9SLeila Ghaffari   CeedVectorDestroy(&qdataSur);
1471*6a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
1472*6a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
1473*6a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
1474*6a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqSur);
1475*6a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxSur);
1476*6a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiSur);
1477*6a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
1478*6a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
1479*6a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
1480*6a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&op_setupSur);
1481*6a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_sur);
1482*6a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_sur);
1483*6a0edaf9SLeila Ghaffari 
1484ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1485ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1486ccaff030SJeremy L Thompson 
1487ccaff030SJeremy L Thompson   // Clean up PETSc
1488ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1489ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1490ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1491ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1492ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1493ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1494ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1495ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1496ccaff030SJeremy L Thompson   return PetscFinalize();
1497ccaff030SJeremy L Thompson }
1498ccaff030SJeremy L Thompson 
1499