xref: /libCEED/examples/solids/src/matops.c (revision 2d93065e886274cb43c7e22f0617419f95f539fb)
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 /// @file
18ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc
19ccaff030SJeremy L Thompson 
20ccaff030SJeremy L Thompson #include "../elasticity.h"
21ccaff030SJeremy L Thompson 
22ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
23ccaff030SJeremy L Thompson // libCEED Operators for MatShell
24ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
25ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
26ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
27ccaff030SJeremy L Thompson   PetscErrorCode ierr;
28ccaff030SJeremy L Thompson   PetscScalar *x, *y;
29ccaff030SJeremy L Thompson 
30ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
31ccaff030SJeremy L Thompson 
32ccaff030SJeremy L Thompson   // Global-to-local
33ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
34ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson   // Setup CEED vectors
37ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr);
38ccaff030SJeremy L Thompson   ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
39ccaff030SJeremy L Thompson   CeedVectorSetArray(user->Xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
40ccaff030SJeremy L Thompson   CeedVectorSetArray(user->Yceed, CEED_MEM_HOST, CEED_USE_POINTER, y);
41ccaff030SJeremy L Thompson 
42ccaff030SJeremy L Thompson   // Apply CEED operator
43ccaff030SJeremy L Thompson   CeedOperatorApply(user->op, user->Xceed, user->Yceed, CEED_REQUEST_IMMEDIATE);
44ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->Yceed, CEED_MEM_HOST);
45ccaff030SJeremy L Thompson 
46ccaff030SJeremy L Thompson   // Restore PETSc vectors
47ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
48ccaff030SJeremy L Thompson   CHKERRQ(ierr);
49ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
50ccaff030SJeremy L Thompson 
51ccaff030SJeremy L Thompson   // Local-to-global
52ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
53ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Yloc, ADD_VALUES, Y); CHKERRQ(ierr);
54ccaff030SJeremy L Thompson 
55ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
56ccaff030SJeremy L Thompson };
57ccaff030SJeremy L Thompson 
58ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
59ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
60ccaff030SJeremy L Thompson   PetscErrorCode ierr;
61ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
62ccaff030SJeremy L Thompson 
63ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
64ccaff030SJeremy L Thompson 
65ccaff030SJeremy L Thompson   // Use computed BCs
66ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
67ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
68ccaff030SJeremy L Thompson                                     user->loadIncrement, NULL, NULL, NULL);
69ccaff030SJeremy L Thompson   CHKERRQ(ierr);
70ccaff030SJeremy L Thompson 
71ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
72ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
73ccaff030SJeremy L Thompson 
74ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
75ccaff030SJeremy L Thompson };
76ccaff030SJeremy L Thompson 
77ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
78ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
79ccaff030SJeremy L Thompson   PetscErrorCode ierr;
80ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
81ccaff030SJeremy L Thompson 
82ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
83ccaff030SJeremy L Thompson 
84ccaff030SJeremy L Thompson   // Use computed BCs
85ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
86ccaff030SJeremy L Thompson 
87ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
88ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
89ccaff030SJeremy L Thompson 
90ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
91ccaff030SJeremy L Thompson };
92ccaff030SJeremy L Thompson 
93ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
94ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
95ccaff030SJeremy L Thompson   PetscErrorCode ierr;
96ccaff030SJeremy L Thompson   UserMult user;
97ccaff030SJeremy L Thompson 
98ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
99ccaff030SJeremy L Thompson 
100ccaff030SJeremy L Thompson   // Zero boundary values
101ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
102ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
103ccaff030SJeremy L Thompson 
104ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
105ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
106ccaff030SJeremy L Thompson 
107ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
108ccaff030SJeremy L Thompson };
109ccaff030SJeremy L Thompson 
110ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
111ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
112ccaff030SJeremy L Thompson   PetscErrorCode ierr;
113ccaff030SJeremy L Thompson   UserMultProlongRestr user;
114ccaff030SJeremy L Thompson   PetscScalar *c, *f;
115ccaff030SJeremy L Thompson 
116ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
117ccaff030SJeremy L Thompson 
118ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
119ccaff030SJeremy L Thompson 
120ccaff030SJeremy L Thompson   // Global-to-local
121ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
122ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmC, X, INSERT_VALUES, user->locVecC);
123ccaff030SJeremy L Thompson   CHKERRQ(ierr);
124ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
125ccaff030SJeremy L Thompson 
126ccaff030SJeremy L Thompson   // Setup CEED vectors
127ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->locVecC, (const PetscScalar **)&c);
128ccaff030SJeremy L Thompson   CHKERRQ(ierr);
129ccaff030SJeremy L Thompson   ierr = VecGetArray(user->locVecF, &f); CHKERRQ(ierr);
130ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, CEED_MEM_HOST, CEED_USE_POINTER, c);
131ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, CEED_MEM_HOST, CEED_USE_POINTER, f);
132ccaff030SJeremy L Thompson 
133ccaff030SJeremy L Thompson   // Apply CEED operator
134ccaff030SJeremy L Thompson   CeedOperatorApply(user->opProlong, user->ceedVecC, user->ceedVecF,
135ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
136ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecF, CEED_MEM_HOST);
137ccaff030SJeremy L Thompson 
138ccaff030SJeremy L Thompson   // Restore PETSc vectors
139ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->locVecC, (const PetscScalar **)c);
140ccaff030SJeremy L Thompson   CHKERRQ(ierr);
141ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->locVecF, &f); CHKERRQ(ierr);
142ccaff030SJeremy L Thompson 
143ccaff030SJeremy L Thompson   // Multiplicity
144ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
145ccaff030SJeremy L Thompson 
146ccaff030SJeremy L Thompson   // Local-to-global
147ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
148ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmF, user->locVecF, ADD_VALUES, Y);
149ccaff030SJeremy L Thompson   CHKERRQ(ierr);
150ccaff030SJeremy L Thompson 
151ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
152ccaff030SJeremy L Thompson }
153ccaff030SJeremy L Thompson 
154ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
155ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
156ccaff030SJeremy L Thompson   PetscErrorCode ierr;
157ccaff030SJeremy L Thompson   UserMultProlongRestr user;
158ccaff030SJeremy L Thompson   PetscScalar *c, *f;
159ccaff030SJeremy L Thompson 
160ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
161ccaff030SJeremy L Thompson 
162ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
163ccaff030SJeremy L Thompson 
164ccaff030SJeremy L Thompson   // Global-to-local
165ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
166ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmF, X, INSERT_VALUES, user->locVecF);
167ccaff030SJeremy L Thompson   CHKERRQ(ierr);
168ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
169ccaff030SJeremy L Thompson 
170ccaff030SJeremy L Thompson   // Multiplicity
171ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
172ccaff030SJeremy L Thompson   CHKERRQ(ierr);
173ccaff030SJeremy L Thompson 
174ccaff030SJeremy L Thompson   // Setup CEED vectors
175ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(user->locVecF, (const PetscScalar **)&f); CHKERRQ(ierr);
176ccaff030SJeremy L Thompson   ierr = VecGetArray(user->locVecC, &c); CHKERRQ(ierr);
177ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, CEED_MEM_HOST, CEED_USE_POINTER, f);
178ccaff030SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, CEED_MEM_HOST, CEED_USE_POINTER, c);
179ccaff030SJeremy L Thompson 
180ccaff030SJeremy L Thompson   // Apply CEED operator
181ccaff030SJeremy L Thompson   CeedOperatorApply(user->opRestrict, user->ceedVecF, user->ceedVecC,
182ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
183ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecC, CEED_MEM_HOST);
184ccaff030SJeremy L Thompson 
185ccaff030SJeremy L Thompson   // Restore PETSc vectors
186ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(user->locVecF, (const PetscScalar **)&f);
187ccaff030SJeremy L Thompson   CHKERRQ(ierr);
188ccaff030SJeremy L Thompson   ierr = VecRestoreArray(user->locVecC, &c); CHKERRQ(ierr);
189ccaff030SJeremy L Thompson 
190ccaff030SJeremy L Thompson   // Local-to-global
191ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
192ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmC, user->locVecC, ADD_VALUES, Y);
193ccaff030SJeremy L Thompson   CHKERRQ(ierr);
194ccaff030SJeremy L Thompson 
195ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
196ccaff030SJeremy L Thompson };
197*2d93065eSjeremylt 
198ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
199ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
200ccaff030SJeremy L Thompson   PetscErrorCode ierr;
201ccaff030SJeremy L Thompson   UserMult user;
202ccaff030SJeremy L Thompson 
203ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
204ccaff030SJeremy L Thompson 
205ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
206ccaff030SJeremy L Thompson 
207ccaff030SJeremy L Thompson   // Compute Diagonal via libCEED
208ccaff030SJeremy L Thompson   CeedVector ceedDiagVec;
209ccaff030SJeremy L Thompson   const CeedScalar *diagArray;
210ccaff030SJeremy L Thompson 
211ccaff030SJeremy L Thompson   // -- Compute Diagonal
212ccaff030SJeremy L Thompson   CeedOperatorAssembleLinearDiagonal(user->op, &ceedDiagVec,
213ccaff030SJeremy L Thompson                                      CEED_REQUEST_IMMEDIATE);
214ccaff030SJeremy L Thompson 
215ccaff030SJeremy L Thompson   // -- Place in PETSc vector
216ccaff030SJeremy L Thompson   CeedVectorGetArrayRead(ceedDiagVec, CEED_MEM_HOST, &diagArray);
217ccaff030SJeremy L Thompson   ierr = VecPlaceArray(user->Xloc, diagArray); CHKERRQ(ierr);
218ccaff030SJeremy L Thompson 
219ccaff030SJeremy L Thompson   // -- Local-to-Global
220ccaff030SJeremy L Thompson   ierr = VecZeroEntries(D); CHKERRQ(ierr);
221ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Xloc, ADD_VALUES, D); CHKERRQ(ierr);
222ccaff030SJeremy L Thompson 
223ccaff030SJeremy L Thompson   // -- Cleanup
224ccaff030SJeremy L Thompson   ierr = VecResetArray(user->Xloc); CHKERRQ(ierr);
225ccaff030SJeremy L Thompson   CeedVectorRestoreArrayRead(ceedDiagVec, &diagArray);
226ccaff030SJeremy L Thompson   CeedVectorDestroy(&ceedDiagVec);
227ccaff030SJeremy L Thompson 
228ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
229ccaff030SJeremy L Thompson };
230*2d93065eSjeremylt 
231*2d93065eSjeremylt // This function calculates the strain energy in the final solution
232*2d93065eSjeremylt PetscErrorCode ComputeStrainEnergy(UserMult user, CeedOperator opEnergy, Vec X,
233*2d93065eSjeremylt                                    CeedVector energyLoc, PetscReal *energy) {
234*2d93065eSjeremylt   PetscErrorCode ierr;
235*2d93065eSjeremylt   PetscScalar *x;
236*2d93065eSjeremylt   CeedInt length;
237*2d93065eSjeremylt 
238*2d93065eSjeremylt   PetscFunctionBeginUser;
239*2d93065eSjeremylt 
240*2d93065eSjeremylt   // Global-to-local
241*2d93065eSjeremylt   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
242*2d93065eSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
243*2d93065eSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
244*2d93065eSjeremylt                                     user->loadIncrement, NULL, NULL, NULL);
245*2d93065eSjeremylt 
246*2d93065eSjeremylt   // Setup libCEED vector
247*2d93065eSjeremylt   ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
248*2d93065eSjeremylt   CHKERRQ(ierr);
249*2d93065eSjeremylt   CeedVectorSetArray(user->Xceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
250*2d93065eSjeremylt 
251*2d93065eSjeremylt   // Apply libCEED operator
252*2d93065eSjeremylt   CeedOperatorApply(opEnergy, user->Xceed, energyLoc, CEED_REQUEST_IMMEDIATE);
253*2d93065eSjeremylt 
254*2d93065eSjeremylt   // Restore PETSc vector
255*2d93065eSjeremylt   ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
256*2d93065eSjeremylt   CHKERRQ(ierr);
257*2d93065eSjeremylt 
258*2d93065eSjeremylt   // Reduce max error
259*2d93065eSjeremylt   const CeedScalar *e;
260*2d93065eSjeremylt   CeedVectorGetArrayRead(energyLoc, CEED_MEM_HOST, &e);
261*2d93065eSjeremylt   (*energy) = 0;
262*2d93065eSjeremylt   CeedVectorGetLength(energyLoc, &length);
263*2d93065eSjeremylt   for (CeedInt i=0; i<length; i++)
264*2d93065eSjeremylt     (*energy) += e[i];
265*2d93065eSjeremylt   CeedVectorRestoreArrayRead(energyLoc, &e);
266*2d93065eSjeremylt 
267*2d93065eSjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
268*2d93065eSjeremylt                        user->comm); CHKERRQ(ierr);
269*2d93065eSjeremylt 
270*2d93065eSjeremylt   PetscFunctionReturn(0);
271*2d93065eSjeremylt };
272