xref: /honee/src/setuplibceed.c (revision 139613f234b672da994bcbdf4b852dbb68ad3ef2)
1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
4a515125bSLeila Ghaffari //
5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and
8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed.
9a515125bSLeila Ghaffari //
10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early
15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari /// @file
18a515125bSLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
19a515125bSLeila Ghaffari 
20a515125bSLeila Ghaffari #include "../navierstokes.h"
21a515125bSLeila Ghaffari 
22a515125bSLeila Ghaffari // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
23a515125bSLeila Ghaffari PetscInt Involute(PetscInt i) {
24a515125bSLeila Ghaffari   return i >= 0 ? i : -(i+1);
25a515125bSLeila Ghaffari }
26a515125bSLeila Ghaffari 
27a515125bSLeila Ghaffari // Utility function to create local CEED restriction
28a515125bSLeila Ghaffari PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
29a515125bSLeila Ghaffari     CeedInt height, DMLabel domain_label,
30a515125bSLeila Ghaffari     CeedInt value, CeedElemRestriction *elem_restr) {
31a515125bSLeila Ghaffari   PetscSection   section;
32a515125bSLeila Ghaffari   PetscInt       p, num_elem, num_dofs, *elem_restrict, elem_offset, num_fields,
33a515125bSLeila Ghaffari                  dim, depth;
34a515125bSLeila Ghaffari   DMLabel        depth_label;
35a515125bSLeila Ghaffari   IS             depth_IS, iter_IS;
36a515125bSLeila Ghaffari   Vec            U_loc;
37a515125bSLeila Ghaffari   const PetscInt *iter_indices;
38a515125bSLeila Ghaffari   PetscErrorCode ierr;
39a515125bSLeila Ghaffari   PetscFunctionBeginUser;
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
42a515125bSLeila Ghaffari   dim -= height;
43a515125bSLeila Ghaffari   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
44a515125bSLeila Ghaffari   ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr);
45a515125bSLeila Ghaffari   PetscInt num_comp[num_fields], field_off[num_fields+1];
46a515125bSLeila Ghaffari   field_off[0] = 0;
47a515125bSLeila Ghaffari   for (PetscInt f=0; f<num_fields; f++) {
48a515125bSLeila Ghaffari     ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr);
49a515125bSLeila Ghaffari     field_off[f+1] = field_off[f] + num_comp[f];
50a515125bSLeila Ghaffari   }
51a515125bSLeila Ghaffari 
52a515125bSLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
53a515125bSLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
54a515125bSLeila Ghaffari   ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_IS);
55a515125bSLeila Ghaffari   CHKERRQ(ierr);
56a515125bSLeila Ghaffari   if (domain_label) {
57a515125bSLeila Ghaffari     IS domain_IS;
58a515125bSLeila Ghaffari     ierr = DMLabelGetStratumIS(domain_label, value, &domain_IS); CHKERRQ(ierr);
59a515125bSLeila Ghaffari     if (domain_IS) { // domain_IS is non-empty
60a515125bSLeila Ghaffari       ierr = ISIntersect(depth_IS, domain_IS, &iter_IS); CHKERRQ(ierr);
61a515125bSLeila Ghaffari       ierr = ISDestroy(&domain_IS); CHKERRQ(ierr);
62a515125bSLeila Ghaffari     } else { // domain_IS is NULL (empty)
63a515125bSLeila Ghaffari       iter_IS = NULL;
64a515125bSLeila Ghaffari     }
65a515125bSLeila Ghaffari     ierr = ISDestroy(&depth_IS); CHKERRQ(ierr);
66a515125bSLeila Ghaffari   } else {
67a515125bSLeila Ghaffari     iter_IS = depth_IS;
68a515125bSLeila Ghaffari   }
69a515125bSLeila Ghaffari   if (iter_IS) {
70a515125bSLeila Ghaffari     ierr = ISGetLocalSize(iter_IS, &num_elem); CHKERRQ(ierr);
71a515125bSLeila Ghaffari     ierr = ISGetIndices(iter_IS, &iter_indices); CHKERRQ(ierr);
72a515125bSLeila Ghaffari   } else {
73a515125bSLeila Ghaffari     num_elem = 0;
74a515125bSLeila Ghaffari     iter_indices = NULL;
75a515125bSLeila Ghaffari   }
76a515125bSLeila Ghaffari   ierr = PetscMalloc1(num_elem*PetscPowInt(P, dim), &elem_restrict);
77a515125bSLeila Ghaffari   CHKERRQ(ierr);
78a515125bSLeila Ghaffari   for (p=0, elem_offset=0; p<num_elem; p++) {
79a515125bSLeila Ghaffari     PetscInt c = iter_indices[p];
80a515125bSLeila Ghaffari     PetscInt num_indices, *indices, num_nodes;
81a515125bSLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
82a515125bSLeila Ghaffari                                    &num_indices, &indices, NULL, NULL);
83a515125bSLeila Ghaffari     CHKERRQ(ierr);
84a515125bSLeila Ghaffari     bool flip = false;
85a515125bSLeila Ghaffari     if (height > 0) {
86a515125bSLeila Ghaffari       PetscInt num_cells, num_faces, start = -1;
87a515125bSLeila Ghaffari       const PetscInt *orients, *faces, *cells;
88a515125bSLeila Ghaffari       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
89a515125bSLeila Ghaffari       ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr);
90a515125bSLeila Ghaffari       if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
91a515125bSLeila Ghaffari                                      "Expected one cell in support of exterior face, but got %D cells",
92a515125bSLeila Ghaffari                                      num_cells);
93a515125bSLeila Ghaffari       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
94a515125bSLeila Ghaffari       ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr);
95a515125bSLeila Ghaffari       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
96a515125bSLeila Ghaffari       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
97a515125bSLeila Ghaffari                                 "Could not find face %D in cone of its support",
98a515125bSLeila Ghaffari                                 c);
99a515125bSLeila Ghaffari       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
100a515125bSLeila Ghaffari       if (orients[start] < 0) flip = true;
101a515125bSLeila Ghaffari     }
102a515125bSLeila Ghaffari     if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF,
103a515125bSLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
104a515125bSLeila Ghaffari           c);
105a515125bSLeila Ghaffari     num_nodes = num_indices / field_off[num_fields];
106a515125bSLeila Ghaffari     for (PetscInt i=0; i<num_nodes; i++) {
107a515125bSLeila Ghaffari       PetscInt ii = i;
108a515125bSLeila Ghaffari       if (flip) {
109a515125bSLeila Ghaffari         if (P == num_nodes) ii = num_nodes - 1 - i;
110a515125bSLeila Ghaffari         else if (P*P == num_nodes) {
111a515125bSLeila Ghaffari           PetscInt row = i / P, col = i % P;
112a515125bSLeila Ghaffari           ii = row + col * P;
113a515125bSLeila Ghaffari         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
114a515125bSLeila Ghaffari                           "No support for flipping point with %D nodes != P (%D) or P^2",
115a515125bSLeila Ghaffari                           num_nodes, P);
116a515125bSLeila Ghaffari       }
117a515125bSLeila Ghaffari       // Check that indices are blocked by node and thus can be coalesced as a single field with
118a515125bSLeila Ghaffari       // field_off[num_fields] = sum(num_comp) components.
119a515125bSLeila Ghaffari       for (PetscInt f=0; f<num_fields; f++) {
120a515125bSLeila Ghaffari         for (PetscInt j=0; j<num_comp[f]; j++) {
121a515125bSLeila Ghaffari           if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j])
122a515125bSLeila Ghaffari               != Involute(indices[ii*num_comp[0]]) + field_off[f] + j)
123a515125bSLeila Ghaffari             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
124a515125bSLeila Ghaffari                      "Cell %D closure indices not interlaced for node %D field %D component %D",
125a515125bSLeila Ghaffari                      c, ii, f, j);
126a515125bSLeila Ghaffari         }
127a515125bSLeila Ghaffari       }
128a515125bSLeila Ghaffari       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
129a515125bSLeila Ghaffari       PetscInt loc = Involute(indices[ii*num_comp[0]]);
130a515125bSLeila Ghaffari       elem_restrict[elem_offset++] = loc;
131a515125bSLeila Ghaffari     }
132a515125bSLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
133a515125bSLeila Ghaffari                                        &num_indices, &indices, NULL, NULL);
134a515125bSLeila Ghaffari     CHKERRQ(ierr);
135a515125bSLeila Ghaffari   }
136a515125bSLeila Ghaffari   if (elem_offset != num_elem*PetscPowInt(P, dim))
137a515125bSLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
138a515125bSLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem,
139a515125bSLeila Ghaffari              PetscPowInt(P, dim),elem_offset);
140a515125bSLeila Ghaffari   if (iter_IS) {
141a515125bSLeila Ghaffari     ierr = ISRestoreIndices(iter_IS, &iter_indices); CHKERRQ(ierr);
142a515125bSLeila Ghaffari   }
143a515125bSLeila Ghaffari   ierr = ISDestroy(&iter_IS); CHKERRQ(ierr);
144a515125bSLeila Ghaffari 
145a515125bSLeila Ghaffari   ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr);
146a515125bSLeila Ghaffari   ierr = VecGetLocalSize(U_loc, &num_dofs); CHKERRQ(ierr);
147a515125bSLeila Ghaffari   ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr);
148a515125bSLeila Ghaffari   CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, dim),
149a515125bSLeila Ghaffari                             field_off[num_fields],
150a515125bSLeila Ghaffari                             1, num_dofs, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restrict,
151a515125bSLeila Ghaffari                             elem_restr);
152a515125bSLeila Ghaffari   ierr = PetscFree(elem_restrict); CHKERRQ(ierr);
153a515125bSLeila Ghaffari   PetscFunctionReturn(0);
154a515125bSLeila Ghaffari }
155a515125bSLeila Ghaffari 
156a515125bSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
157a515125bSLeila Ghaffari PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
158a515125bSLeila Ghaffari                                        DMLabel domain_label, PetscInt value,
159a515125bSLeila Ghaffari                                        CeedInt P, CeedInt Q, CeedInt q_data_size,
160a515125bSLeila Ghaffari                                        CeedElemRestriction *elem_restr_q,
161a515125bSLeila Ghaffari                                        CeedElemRestriction *elem_restr_x,
162a515125bSLeila Ghaffari                                        CeedElemRestriction *elem_restr_qd_i) {
163a515125bSLeila Ghaffari   DM             dm_coord;
164a515125bSLeila Ghaffari   CeedInt        dim, loc_num_elem;
165a515125bSLeila Ghaffari   CeedInt        Q_dim;
166a515125bSLeila Ghaffari   PetscErrorCode ierr;
167a515125bSLeila Ghaffari   PetscFunctionBeginUser;
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
170a515125bSLeila Ghaffari   dim -= height;
171a515125bSLeila Ghaffari   Q_dim = CeedIntPow(Q, dim);
172a515125bSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
173a515125bSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
174a515125bSLeila Ghaffari   CHKERRQ(ierr);
175a515125bSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domain_label, value,
176*139613f2SLeila Ghaffari                                    elem_restr_q); CHKERRQ(ierr);
177a515125bSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm_coord, 2, height, domain_label, value,
178*139613f2SLeila Ghaffari                                    elem_restr_x); CHKERRQ(ierr);
179a515125bSLeila Ghaffari   CeedElemRestrictionGetNumElements(*elem_restr_q, &loc_num_elem);
180a515125bSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q_dim,
181a515125bSLeila Ghaffari                                    q_data_size, q_data_size*loc_num_elem*Q_dim,
182a515125bSLeila Ghaffari                                    CEED_STRIDES_BACKEND, elem_restr_qd_i);
183a515125bSLeila Ghaffari   PetscFunctionReturn(0);
184a515125bSLeila Ghaffari }
185a515125bSLeila Ghaffari 
186a515125bSLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
187a515125bSLeila Ghaffari PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
188a515125bSLeila Ghaffari                                        CeedData ceed_data, Physics phys,
189a515125bSLeila Ghaffari                                        CeedOperator op_apply_vol, CeedInt height,
190a515125bSLeila Ghaffari                                        CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
191a515125bSLeila Ghaffari                                        CeedOperator *op_apply) {
192a515125bSLeila Ghaffari   CeedInt        dim, num_face;
193a515125bSLeila Ghaffari   DMLabel        domain_label;
194a515125bSLeila Ghaffari   PetscErrorCode ierr;
195a515125bSLeila Ghaffari   PetscFunctionBeginUser;
196a515125bSLeila Ghaffari 
197a515125bSLeila Ghaffari   // Create Composite Operaters
198a515125bSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
199a515125bSLeila Ghaffari 
200a515125bSLeila Ghaffari   // --Apply Sub-Operator for the volume
201a515125bSLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_apply_vol);
202a515125bSLeila Ghaffari 
203a515125bSLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
204*139613f2SLeila Ghaffari   if (phys->has_neumann) {
205a515125bSLeila Ghaffari     // --- Setup
206a515125bSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domain_label); CHKERRQ(ierr);
207a515125bSLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
208a515125bSLeila Ghaffari     if (dim == 2) num_face = 4;
209a515125bSLeila Ghaffari     if (dim == 3) num_face = 6;
210a515125bSLeila Ghaffari 
211a515125bSLeila Ghaffari     // --- Get number of quadrature points for the boundaries
212a515125bSLeila Ghaffari     CeedInt num_qpts_sur;
213a515125bSLeila Ghaffari     CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur);
214a515125bSLeila Ghaffari 
215a515125bSLeila Ghaffari     // --- Create Sub-Operator for each face
216a515125bSLeila Ghaffari     for (CeedInt i=0; i<num_face; i++) {
217a515125bSLeila Ghaffari       CeedVector          q_data_sur;
218a515125bSLeila Ghaffari       CeedOperator        op_setup_sur, op_apply_sur;
219a515125bSLeila Ghaffari       CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur,
220a515125bSLeila Ghaffari                           elem_restr_qd_i_sur;
221a515125bSLeila Ghaffari 
222a515125bSLeila Ghaffari       // ---- CEED Restriction
223a515125bSLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domain_label, i+1, P_sur,
224a515125bSLeila Ghaffari                                      Q_sur, q_data_size_sur, &elem_restr_q_sur,
225a515125bSLeila Ghaffari                                      &elem_restr_x_sur, &elem_restr_qd_i_sur);
226a515125bSLeila Ghaffari       CHKERRQ(ierr);
227a515125bSLeila Ghaffari 
228a515125bSLeila Ghaffari       // ---- CEED Vector
229a515125bSLeila Ghaffari       PetscInt loc_num_elem_sur;
230a515125bSLeila Ghaffari       CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur);
231a515125bSLeila Ghaffari       CeedVectorCreate(ceed, q_data_size_sur*loc_num_elem_sur*num_qpts_sur,
232a515125bSLeila Ghaffari                        &q_data_sur);
233a515125bSLeila Ghaffari 
234a515125bSLeila Ghaffari       // ---- CEED Operator
235a515125bSLeila Ghaffari       // ----- CEED Operator for Setup (geometric factors)
236a515125bSLeila Ghaffari       CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
237a515125bSLeila Ghaffari       CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur,
238*139613f2SLeila Ghaffari                            ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE);
239a515125bSLeila Ghaffari       CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE,
240a515125bSLeila Ghaffari                            ceed_data->basis_x_sur, CEED_VECTOR_NONE);
241a515125bSLeila Ghaffari       CeedOperatorSetField(op_setup_sur, "q_data_sur", elem_restr_qd_i_sur,
242a515125bSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
243a515125bSLeila Ghaffari 
244a515125bSLeila Ghaffari       // ----- CEED Operator for Physics
245a515125bSLeila Ghaffari       CeedOperatorCreate(ceed, ceed_data->qf_apply_sur, NULL, NULL, &op_apply_sur);
246a515125bSLeila Ghaffari       CeedOperatorSetField(op_apply_sur, "q", elem_restr_q_sur,
247a515125bSLeila Ghaffari                            ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
248a515125bSLeila Ghaffari       CeedOperatorSetField(op_apply_sur, "q_data_sur", elem_restr_qd_i_sur,
249a515125bSLeila Ghaffari                            CEED_BASIS_COLLOCATED, q_data_sur);
250a515125bSLeila Ghaffari       CeedOperatorSetField(op_apply_sur, "x", elem_restr_x_sur,
251a515125bSLeila Ghaffari                            ceed_data->basis_x_sur, ceed_data->x_coord);
252a515125bSLeila Ghaffari       CeedOperatorSetField(op_apply_sur, "v", elem_restr_q_sur,
253a515125bSLeila Ghaffari                            ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
254a515125bSLeila Ghaffari 
255a515125bSLeila Ghaffari       // ----- Apply CEED operator for Setup
256a515125bSLeila Ghaffari       CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur,
257a515125bSLeila Ghaffari                         CEED_REQUEST_IMMEDIATE);
258a515125bSLeila Ghaffari 
259a515125bSLeila Ghaffari       // ----- Apply Sub-Operator for the Boundary
260a515125bSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_apply_sur);
261a515125bSLeila Ghaffari 
262a515125bSLeila Ghaffari       // ----- Cleanup
263a515125bSLeila Ghaffari       CeedVectorDestroy(&q_data_sur);
264a515125bSLeila Ghaffari       CeedElemRestrictionDestroy(&elem_restr_q_sur);
265a515125bSLeila Ghaffari       CeedElemRestrictionDestroy(&elem_restr_x_sur);
266a515125bSLeila Ghaffari       CeedElemRestrictionDestroy(&elem_restr_qd_i_sur);
267a515125bSLeila Ghaffari       CeedOperatorDestroy(&op_setup_sur);
268a515125bSLeila Ghaffari       CeedOperatorDestroy(&op_apply_sur);
269a515125bSLeila Ghaffari     }
270a515125bSLeila Ghaffari   }
271a515125bSLeila Ghaffari 
272a515125bSLeila Ghaffari   PetscFunctionReturn(0);
273a515125bSLeila Ghaffari }
274a515125bSLeila Ghaffari 
275a515125bSLeila Ghaffari PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user,
276a515125bSLeila Ghaffari                             AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
277a515125bSLeila Ghaffari   PetscErrorCode ierr;
278a515125bSLeila Ghaffari   PetscFunctionBeginUser;
279a515125bSLeila Ghaffari 
280a515125bSLeila Ghaffari   // *****************************************************************************
281a515125bSLeila Ghaffari   // Set up CEED objects for the interior domain (volume)
282a515125bSLeila Ghaffari   // *****************************************************************************
283a515125bSLeila Ghaffari   const PetscInt num_comp_q      = 5;
284a515125bSLeila Ghaffari   const CeedInt  dim             = problem->dim,
285a515125bSLeila Ghaffari                  num_comp_x      = problem->dim,
286a515125bSLeila Ghaffari                  q_data_size_vol = problem->q_data_size_vol,
287a515125bSLeila Ghaffari                  P               = app_ctx->degree + 1,
288a515125bSLeila Ghaffari                  Q               = P + app_ctx->q_extra;
289a515125bSLeila Ghaffari 
290a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
291a515125bSLeila Ghaffari   // CEED Bases
292a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
293a515125bSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS,
294a515125bSLeila Ghaffari                                   &ceed_data->basis_q);
295a515125bSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS,
296a515125bSLeila Ghaffari                                   &ceed_data->basis_x);
297a515125bSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P,
298a515125bSLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &ceed_data->basis_xc);
299a515125bSLeila Ghaffari 
300a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
301a515125bSLeila Ghaffari   // CEED Restrictions
302a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
303a515125bSLeila Ghaffari   // -- Create restriction
304a515125bSLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, P, Q,
305a515125bSLeila Ghaffari                                  q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x,
306a515125bSLeila Ghaffari                                  &ceed_data->elem_restr_qd_i); CHKERRQ(ierr);
307a515125bSLeila Ghaffari   // -- Create E vectors
308a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL);
309a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed,
310a515125bSLeila Ghaffari                                   NULL);
311a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL);
312a515125bSLeila Ghaffari 
313a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
314a515125bSLeila Ghaffari   // CEED QFunctions
315a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
316a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
317a515125bSLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol, problem->setup_vol_loc,
318a515125bSLeila Ghaffari                               &ceed_data->qf_setup_vol);
319a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x*dim,
320a515125bSLeila Ghaffari                         CEED_EVAL_GRAD);
321a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT);
322a515125bSLeila Ghaffari   CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "q_data", q_data_size_vol,
323a515125bSLeila Ghaffari                          CEED_EVAL_NONE);
324a515125bSLeila Ghaffari 
325a515125bSLeila Ghaffari   // -- Create QFunction for ICs
326a515125bSLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc,
327a515125bSLeila Ghaffari                               &ceed_data->qf_ics);
328a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP);
329a515125bSLeila Ghaffari   CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE);
330a515125bSLeila Ghaffari 
331a515125bSLeila Ghaffari   // -- Create QFunction for RHS
332a515125bSLeila Ghaffari   if (problem->apply_vol_rhs) {
333a515125bSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs,
334a515125bSLeila Ghaffari                                 problem->apply_vol_rhs_loc, &ceed_data->qf_rhs_vol);
335a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP);
336a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "dq", num_comp_q*dim,
337a515125bSLeila Ghaffari                           CEED_EVAL_GRAD);
338a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q_data", q_data_size_vol,
339a515125bSLeila Ghaffari                           CEED_EVAL_NONE);
340a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP);
341a515125bSLeila Ghaffari     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q,
342a515125bSLeila Ghaffari                            CEED_EVAL_INTERP);
343a515125bSLeila Ghaffari     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "dv", num_comp_q*dim,
344a515125bSLeila Ghaffari                            CEED_EVAL_GRAD);
345a515125bSLeila Ghaffari   }
346a515125bSLeila Ghaffari 
347a515125bSLeila Ghaffari   // -- Create QFunction for IFunction
348a515125bSLeila Ghaffari   if (problem->apply_vol_ifunction) {
349a515125bSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction,
350a515125bSLeila Ghaffari                                 problem->apply_vol_ifunction_loc, &ceed_data->qf_ifunction_vol);
351a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q,
352a515125bSLeila Ghaffari                           CEED_EVAL_INTERP);
353a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "dq", num_comp_q*dim,
354a515125bSLeila Ghaffari                           CEED_EVAL_GRAD);
355a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q_dot", num_comp_q,
356a515125bSLeila Ghaffari                           CEED_EVAL_INTERP);
357a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q_data", q_data_size_vol,
358a515125bSLeila Ghaffari                           CEED_EVAL_NONE);
359a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x,
360a515125bSLeila Ghaffari                           CEED_EVAL_INTERP);
361a515125bSLeila Ghaffari     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q,
362a515125bSLeila Ghaffari                            CEED_EVAL_INTERP);
363a515125bSLeila Ghaffari     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "dv", num_comp_q*dim,
364a515125bSLeila Ghaffari                            CEED_EVAL_GRAD);
365a515125bSLeila Ghaffari   }
366a515125bSLeila Ghaffari 
367a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
368a515125bSLeila Ghaffari   // Element coordinates
369a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
370a515125bSLeila Ghaffari   // -- Create CEED vector
371a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord,
372a515125bSLeila Ghaffari                                   NULL);
373a515125bSLeila Ghaffari 
374a515125bSLeila Ghaffari   // -- Copy PETSc vector in CEED vector
375a515125bSLeila Ghaffari   Vec               X_loc;
376a515125bSLeila Ghaffari   const PetscScalar *X_loc_array;
377a515125bSLeila Ghaffari   ierr = DMGetCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);
378a515125bSLeila Ghaffari   ierr = VecGetArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr);
379a515125bSLeila Ghaffari   CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES,
380a515125bSLeila Ghaffari                      (PetscScalar *)X_loc_array);
381a515125bSLeila Ghaffari   ierr = VecRestoreArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr);
382a515125bSLeila Ghaffari 
383a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
384a515125bSLeila Ghaffari   // CEED vectors
385a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
386a515125bSLeila Ghaffari   // -- Create CEED vector for geometric data
387a515125bSLeila Ghaffari   CeedInt  num_qpts_vol;
388a515125bSLeila Ghaffari   PetscInt loc_num_elem_vol;
389a515125bSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts_vol);
390a515125bSLeila Ghaffari   CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &loc_num_elem_vol);
391a515125bSLeila Ghaffari   CeedVectorCreate(ceed, q_data_size_vol*loc_num_elem_vol*num_qpts_vol,
392a515125bSLeila Ghaffari                    &ceed_data->q_data);
393a515125bSLeila Ghaffari 
394a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
395a515125bSLeila Ghaffari   // CEED Operators
396a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
397a515125bSLeila Ghaffari   // -- Create CEED operator for quadrature data
398a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL,
399a515125bSLeila Ghaffari                      &ceed_data->op_setup_vol);
400a515125bSLeila Ghaffari   CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x,
401a515125bSLeila Ghaffari                        ceed_data->basis_x, CEED_VECTOR_ACTIVE);
402a515125bSLeila Ghaffari   CeedOperatorSetField(ceed_data->op_setup_vol, "weight",
403*139613f2SLeila Ghaffari                        CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE);
404a515125bSLeila Ghaffari   CeedOperatorSetField(ceed_data->op_setup_vol, "q_data",
405*139613f2SLeila Ghaffari                        ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
406a515125bSLeila Ghaffari 
407a515125bSLeila Ghaffari   // -- Create CEED operator for ICs
408a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics);
409a515125bSLeila Ghaffari   CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x,
410a515125bSLeila Ghaffari                        ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
411a515125bSLeila Ghaffari   CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q,
412a515125bSLeila Ghaffari                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
413a515125bSLeila Ghaffari 
414a515125bSLeila Ghaffari   // Create CEED operator for RHS
415a515125bSLeila Ghaffari   if (ceed_data->qf_rhs_vol) {
416a515125bSLeila Ghaffari     CeedOperator op;
417a515125bSLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op);
418a515125bSLeila Ghaffari     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
419a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
420a515125bSLeila Ghaffari     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q,
421a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
422a515125bSLeila Ghaffari     CeedOperatorSetField(op, "q_data", ceed_data->elem_restr_qd_i,
423*139613f2SLeila Ghaffari                          CEED_BASIS_COLLOCATED, ceed_data->q_data);
424a515125bSLeila Ghaffari     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x,
425a515125bSLeila Ghaffari                          ceed_data->x_coord);
426a515125bSLeila Ghaffari     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
427a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
428a515125bSLeila Ghaffari     CeedOperatorSetField(op, "dv", ceed_data->elem_restr_q, ceed_data->basis_q,
429a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
430a515125bSLeila Ghaffari     user->op_rhs_vol = op;
431a515125bSLeila Ghaffari   }
432a515125bSLeila Ghaffari 
433a515125bSLeila Ghaffari   // -- CEED operator for IFunction
434a515125bSLeila Ghaffari   if (ceed_data->qf_ifunction_vol) {
435a515125bSLeila Ghaffari     CeedOperator op;
436a515125bSLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op);
437a515125bSLeila Ghaffari     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
438a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
439a515125bSLeila Ghaffari     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q,
440a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
441a515125bSLeila Ghaffari     CeedOperatorSetField(op, "q_dot", ceed_data->elem_restr_q, ceed_data->basis_q,
442a515125bSLeila Ghaffari                          user->q_dot_ceed);
443a515125bSLeila Ghaffari     CeedOperatorSetField(op, "q_data", ceed_data->elem_restr_qd_i,
444*139613f2SLeila Ghaffari                          CEED_BASIS_COLLOCATED, ceed_data->q_data);
445a515125bSLeila Ghaffari     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x,
446a515125bSLeila Ghaffari                          ceed_data->x_coord);
447a515125bSLeila Ghaffari     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
448a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
449a515125bSLeila Ghaffari     CeedOperatorSetField(op, "dv", ceed_data->elem_restr_q, ceed_data->basis_q,
450a515125bSLeila Ghaffari                          CEED_VECTOR_ACTIVE);
451a515125bSLeila Ghaffari     user->op_ifunction_vol = op;
452a515125bSLeila Ghaffari   }
453a515125bSLeila Ghaffari 
454a515125bSLeila Ghaffari   // *****************************************************************************
455a515125bSLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
456a515125bSLeila Ghaffari   // *****************************************************************************
457a515125bSLeila Ghaffari   CeedInt height  = 1,
458a515125bSLeila Ghaffari           dim_sur = dim - height,
459a515125bSLeila Ghaffari           P_sur   = app_ctx->degree + 1,
460a515125bSLeila Ghaffari           Q_sur   = P_sur + app_ctx->q_extra;
461a515125bSLeila Ghaffari   const CeedInt q_data_size_sur = problem->q_data_size_sur;
462a515125bSLeila Ghaffari 
463a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
464a515125bSLeila Ghaffari   // CEED Bases
465a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
466a515125bSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur,
467a515125bSLeila Ghaffari                                   CEED_GAUSS, &ceed_data->basis_q_sur);
468a515125bSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS,
469a515125bSLeila Ghaffari                                   &ceed_data->basis_x_sur);
470a515125bSLeila Ghaffari 
471a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
472a515125bSLeila Ghaffari   // CEED QFunctions
473a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
474a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
475a515125bSLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur, problem->setup_sur_loc,
476a515125bSLeila Ghaffari                               &ceed_data->qf_setup_sur);
477a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x*dim_sur,
478a515125bSLeila Ghaffari                         CEED_EVAL_GRAD);
479a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT);
480a515125bSLeila Ghaffari   CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "q_data_sur", q_data_size_sur,
481a515125bSLeila Ghaffari                          CEED_EVAL_NONE);
482a515125bSLeila Ghaffari 
483a515125bSLeila Ghaffari   // -- Creat QFunction for the physics on the boundaries
484a515125bSLeila Ghaffari   if (problem->apply_sur) {
485a515125bSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->apply_sur, problem->apply_sur_loc,
486a515125bSLeila Ghaffari                                 &ceed_data->qf_apply_sur);
487a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_apply_sur, "q", num_comp_q,
488a515125bSLeila Ghaffari                           CEED_EVAL_INTERP);
489a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_apply_sur, "q_data_sur", q_data_size_sur,
490a515125bSLeila Ghaffari                           CEED_EVAL_NONE);
491a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_apply_sur, "x", num_comp_x,
492a515125bSLeila Ghaffari                           CEED_EVAL_INTERP);
493a515125bSLeila Ghaffari     CeedQFunctionAddOutput(ceed_data->qf_apply_sur, "v", num_comp_q,
494a515125bSLeila Ghaffari                            CEED_EVAL_INTERP);
495a515125bSLeila Ghaffari   }
496a515125bSLeila Ghaffari 
497a515125bSLeila Ghaffari   // *****************************************************************************
498a515125bSLeila Ghaffari   // CEED Operator Apply
499a515125bSLeila Ghaffari   // *****************************************************************************
500a515125bSLeila Ghaffari   // -- Apply CEED Operator for the geometric data
501a515125bSLeila Ghaffari   CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord,
502a515125bSLeila Ghaffari                     ceed_data->q_data, CEED_REQUEST_IMMEDIATE);
503a515125bSLeila Ghaffari 
504a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
505a515125bSLeila Ghaffari   if (!user->phys->implicit) { // RHS
506a515125bSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys,
507a515125bSLeila Ghaffari                                    user->op_rhs_vol, height, P_sur, Q_sur,
508a515125bSLeila Ghaffari                                    q_data_size_sur, &user->op_rhs); CHKERRQ(ierr);
509a515125bSLeila Ghaffari   } else { // IFunction
510a515125bSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys,
511a515125bSLeila Ghaffari                                    user->op_ifunction_vol, height, P_sur, Q_sur,
512a515125bSLeila Ghaffari                                    q_data_size_sur, &user->op_ifunction); CHKERRQ(ierr);
513a515125bSLeila Ghaffari   }
514a515125bSLeila Ghaffari 
515a515125bSLeila Ghaffari   PetscFunctionReturn(0);
516a515125bSLeila Ghaffari }
517