xref: /libCEED/examples/solids/src/matops.c (revision 2fcdcaf2bb91d66e7e43fbde06edc0ff61e53133)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 /// @file
18 /// Matrix shell operations for solid mechanics example using PETSc
19 
20 #include "../elasticity.h"
21 
22 // -----------------------------------------------------------------------------
23 // libCEED Operators for MatShell
24 // -----------------------------------------------------------------------------
25 // This function uses libCEED to compute the local action of an operator
26 PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
27   PetscErrorCode ierr;
28   PetscScalar *x, *y;
29 
30   PetscFunctionBeginUser;
31 
32   // Global-to-local
33   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
34   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
35 
36   // Setup CEED vectors
37   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
38   CHKERRQ(ierr);
39   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
40   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
41   CeedVectorSetArray(user->Yceed, user->memType, CEED_USE_POINTER, y);
42 
43   // Apply CEED operator
44   // Note: We could use VecGetArrayInPlace. Instead, we use SetArray/TakeArray
45   //         so we can request host memory for easier debugging.
46   CeedOperatorApply(user->op, user->Xceed, user->Yceed, CEED_REQUEST_IMMEDIATE);
47 
48   // Restore PETSc vectors
49   CeedVectorTakeArray(user->Xceed, user->memType, NULL);
50   CeedVectorTakeArray(user->Yceed, user->memType, NULL);
51   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
52   CHKERRQ(ierr);
53   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
54 
55   // Local-to-global
56   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
57   ierr = DMLocalToGlobal(user->dm, user->Yloc, ADD_VALUES, Y); CHKERRQ(ierr);
58 
59   PetscFunctionReturn(0);
60 };
61 
62 // This function uses libCEED to compute the non-linear residual
63 PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
64   PetscErrorCode ierr;
65   UserMult user = (UserMult)ctx;
66 
67   PetscFunctionBeginUser;
68 
69   // Use computed BCs
70   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
71   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
72                                     user->loadIncrement, NULL, NULL, NULL);
73   CHKERRQ(ierr);
74 
75   // libCEED for local action of residual evaluator
76   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
77 
78   PetscFunctionReturn(0);
79 };
80 
81 // This function uses libCEED to apply the Jacobian for assembly via a SNES
82 PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
83   PetscErrorCode ierr;
84   UserMult user = (UserMult)ctx;
85 
86   PetscFunctionBeginUser;
87 
88   // Zero boundary values
89   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
90 
91   // libCEED for local action of residual evaluator
92   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
93 
94   PetscFunctionReturn(0);
95 };
96 
97 // This function uses libCEED to compute the action of the Jacobian
98 PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
99   PetscErrorCode ierr;
100   UserMult user;
101 
102   PetscFunctionBeginUser;
103 
104   // Zero boundary values
105   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
106   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
107 
108   // libCEED for local action of Jacobian
109   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
110 
111   PetscFunctionReturn(0);
112 };
113 
114 // This function uses libCEED to compute the action of the prolongation operator
115 PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
116   PetscErrorCode ierr;
117   UserMultProlongRestr user;
118   PetscScalar *c, *f;
119 
120   PetscFunctionBeginUser;
121 
122   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
123 
124   // Global-to-local
125   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
126   ierr = DMGlobalToLocal(user->dmC, X, INSERT_VALUES, user->locVecC);
127   CHKERRQ(ierr);
128   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
129 
130   // Setup CEED vectors
131   ierr = user->VecGetArrayRead(user->locVecC, (const PetscScalar **)&c);
132   CHKERRQ(ierr);
133   ierr = user->VecGetArray(user->locVecF, &f); CHKERRQ(ierr);
134   CeedVectorSetArray(user->ceedVecC, user->memType, CEED_USE_POINTER, c);
135   CeedVectorSetArray(user->ceedVecF, user->memType, CEED_USE_POINTER, f);
136 
137   // Apply CEED operator
138   CeedOperatorApply(user->opProlong, user->ceedVecC, user->ceedVecF,
139                     CEED_REQUEST_IMMEDIATE);
140 
141   // Restore PETSc vectors
142   CeedVectorTakeArray(user->ceedVecC, user->memType, NULL);
143   CeedVectorTakeArray(user->ceedVecF, user->memType, NULL);
144   ierr = user->VecRestoreArrayRead(user->locVecC, (const PetscScalar **)&c);
145   CHKERRQ(ierr);
146   ierr = user->VecRestoreArray(user->locVecF, &f); CHKERRQ(ierr);
147 
148   // Multiplicity
149   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
150 
151   // Local-to-global
152   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
153   ierr = DMLocalToGlobal(user->dmF, user->locVecF, ADD_VALUES, Y);
154   CHKERRQ(ierr);
155 
156   PetscFunctionReturn(0);
157 }
158 
159 // This function uses libCEED to compute the action of the restriction operator
160 PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
161   PetscErrorCode ierr;
162   UserMultProlongRestr user;
163   PetscScalar *c, *f;
164 
165   PetscFunctionBeginUser;
166 
167   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
168 
169   // Global-to-local
170   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
171   ierr = DMGlobalToLocal(user->dmF, X, INSERT_VALUES, user->locVecF);
172   CHKERRQ(ierr);
173   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
174 
175   // Multiplicity
176   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
177   CHKERRQ(ierr);
178 
179   // Setup CEED vectors
180   ierr = user->VecGetArrayRead(user->locVecF, (const PetscScalar **)&f);
181   CHKERRQ(ierr);
182   ierr = user->VecGetArray(user->locVecC, &c); CHKERRQ(ierr);
183   CeedVectorSetArray(user->ceedVecF, user->memType, CEED_USE_POINTER, f);
184   CeedVectorSetArray(user->ceedVecC, user->memType, CEED_USE_POINTER, c);
185 
186   // Apply CEED operator
187   CeedOperatorApply(user->opRestrict, user->ceedVecF, user->ceedVecC,
188                     CEED_REQUEST_IMMEDIATE);
189 
190   // Restore PETSc vectors
191   CeedVectorTakeArray(user->ceedVecF, user->memType, NULL);
192   CeedVectorTakeArray(user->ceedVecC, user->memType, NULL);
193   ierr = user->VecRestoreArrayRead(user->locVecF, (const PetscScalar **)&f);
194   CHKERRQ(ierr);
195   ierr = user->VecRestoreArray(user->locVecC, &c); CHKERRQ(ierr);
196 
197   // Local-to-global
198   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
199   ierr = DMLocalToGlobal(user->dmC, user->locVecC, ADD_VALUES, Y);
200   CHKERRQ(ierr);
201 
202   PetscFunctionReturn(0);
203 };
204 
205 // This function returns the computed diagonal of the operator
206 PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
207   PetscErrorCode ierr;
208   UserMult user;
209 
210   PetscFunctionBeginUser;
211 
212   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
213 
214   // -- Set physics context
215   if (user->physSmoother)
216     CeedQFunctionSetContext(user->qf, user->physSmoother,
217                             sizeof(*user->physSmoother));
218 
219   // Compute Diagonal via libCEED
220   PetscScalar *x;
221 
222   // -- Place PETSc vector in libCEED vector
223   ierr = user->VecGetArray(user->Xloc, &x); CHKERRQ(ierr);
224   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
225 
226   // -- Compute Diagonal
227   CeedOperatorLinearAssembleDiagonal(user->op, user->Xceed,
228                                      CEED_REQUEST_IMMEDIATE);
229 
230   // -- Reset physics context
231   if (user->physSmoother)
232     CeedQFunctionSetContext(user->qf, user->phys, sizeof(*user->phys));
233 
234   // -- Local-to-Global
235   CeedVectorTakeArray(user->Xceed, user->memType, NULL);
236   ierr = user->VecRestoreArray(user->Xloc, &x); CHKERRQ(ierr);
237   ierr = VecZeroEntries(D); CHKERRQ(ierr);
238   ierr = DMLocalToGlobal(user->dm, user->Xloc, ADD_VALUES, D); CHKERRQ(ierr);
239 
240   // Cleanup
241   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
242 
243   PetscFunctionReturn(0);
244 };
245 
246 // This function calculates the strain energy in the final solution
247 PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
248                                    CeedOperator opEnergy, Vec X,
249                                    PetscReal *energy) {
250   PetscErrorCode ierr;
251   PetscScalar *x;
252   CeedInt length;
253 
254   PetscFunctionBeginUser;
255 
256   // Global-to-local
257   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
258   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
259   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
260                                     user->loadIncrement, NULL, NULL, NULL);
261   CHKERRQ(ierr);
262 
263   // Setup libCEED input vector
264   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
265   CHKERRQ(ierr);
266   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
267 
268   // Setup libCEED output vector
269   Vec Eloc;
270   CeedVector eloc;
271   ierr = DMCreateLocalVector(dmEnergy, &Eloc); CHKERRQ(ierr);
272   ierr = VecGetSize(Eloc, &length); CHKERRQ(ierr);
273   ierr = VecDestroy(&Eloc); CHKERRQ(ierr);
274   CeedVectorCreate(user->ceed, length, &eloc);
275 
276   // Apply libCEED operator
277   CeedOperatorApply(opEnergy, user->Xceed, eloc, CEED_REQUEST_IMMEDIATE);
278 
279   // Restore PETSc vector
280   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
281   CHKERRQ(ierr);
282 
283   // Reduce max error
284   const CeedScalar *e;
285   CeedVectorGetArrayRead(eloc, CEED_MEM_HOST, &e);
286   (*energy) = 0;
287   for (CeedInt i=0; i<length; i++)
288     (*energy) += e[i];
289   CeedVectorRestoreArrayRead(eloc, &e);
290   CeedVectorDestroy(&eloc);
291 
292   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
293                        user->comm); CHKERRQ(ierr);
294 
295   PetscFunctionReturn(0);
296 };
297