xref: /libCEED/examples/solids/src/matops.c (revision ee83402c013a7f981c352f7f4774f539ff3e4da7)
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   // Local-to-global
149   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
150   ierr = DMLocalToGlobal(user->dmF, user->locVecF, ADD_VALUES, Y);
151   CHKERRQ(ierr);
152 
153   PetscFunctionReturn(0);
154 }
155 
156 // This function uses libCEED to compute the action of the restriction operator
157 PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
158   PetscErrorCode ierr;
159   UserMultProlongRestr user;
160   PetscScalar *c, *f;
161 
162   PetscFunctionBeginUser;
163 
164   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
165 
166   // Global-to-local
167   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
168   ierr = DMGlobalToLocal(user->dmF, X, INSERT_VALUES, user->locVecF);
169   CHKERRQ(ierr);
170   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
171 
172   // Setup CEED vectors
173   ierr = user->VecGetArrayRead(user->locVecF, (const PetscScalar **)&f);
174   CHKERRQ(ierr);
175   ierr = user->VecGetArray(user->locVecC, &c); CHKERRQ(ierr);
176   CeedVectorSetArray(user->ceedVecF, user->memType, CEED_USE_POINTER, f);
177   CeedVectorSetArray(user->ceedVecC, user->memType, CEED_USE_POINTER, c);
178 
179   // Apply CEED operator
180   CeedOperatorApply(user->opRestrict, user->ceedVecF, user->ceedVecC,
181                     CEED_REQUEST_IMMEDIATE);
182 
183   // Restore PETSc vectors
184   CeedVectorTakeArray(user->ceedVecF, user->memType, NULL);
185   CeedVectorTakeArray(user->ceedVecC, user->memType, NULL);
186   ierr = user->VecRestoreArrayRead(user->locVecF, (const PetscScalar **)&f);
187   CHKERRQ(ierr);
188   ierr = user->VecRestoreArray(user->locVecC, &c); CHKERRQ(ierr);
189 
190   // Local-to-global
191   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
192   ierr = DMLocalToGlobal(user->dmC, user->locVecC, ADD_VALUES, Y);
193   CHKERRQ(ierr);
194 
195   PetscFunctionReturn(0);
196 };
197 
198 // This function returns the computed diagonal of the operator
199 PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
200   PetscErrorCode ierr;
201   UserMult user;
202 
203   PetscFunctionBeginUser;
204 
205   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
206 
207   // -- Set physics context
208   if (user->ctxPhysSmoother)
209     CeedQFunctionSetContext(user->qf, user->ctxPhysSmoother);
210 
211   // Compute Diagonal via libCEED
212   PetscScalar *x;
213 
214   // -- Place PETSc vector in libCEED vector
215   ierr = user->VecGetArray(user->Xloc, &x); CHKERRQ(ierr);
216   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
217 
218   // -- Compute Diagonal
219   CeedOperatorLinearAssembleDiagonal(user->op, user->Xceed,
220                                      CEED_REQUEST_IMMEDIATE);
221 
222   // -- Reset physics context
223   if (user->ctxPhysSmoother)
224     CeedQFunctionSetContext(user->qf, user->ctxPhys);
225 
226   // -- Local-to-Global
227   CeedVectorTakeArray(user->Xceed, user->memType, NULL);
228   ierr = user->VecRestoreArray(user->Xloc, &x); CHKERRQ(ierr);
229   ierr = VecZeroEntries(D); CHKERRQ(ierr);
230   ierr = DMLocalToGlobal(user->dm, user->Xloc, ADD_VALUES, D); CHKERRQ(ierr);
231 
232   // Cleanup
233   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
234 
235   PetscFunctionReturn(0);
236 };
237 
238 // This function calculates the strain energy in the final solution
239 PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
240                                    CeedOperator opEnergy, Vec X,
241                                    PetscReal *energy) {
242   PetscErrorCode ierr;
243   PetscScalar *x;
244   CeedInt length;
245 
246   PetscFunctionBeginUser;
247 
248   // Global-to-local
249   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
250   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
251   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
252                                     user->loadIncrement, NULL, NULL, NULL);
253   CHKERRQ(ierr);
254 
255   // Setup libCEED input vector
256   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
257   CHKERRQ(ierr);
258   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
259 
260   // Setup libCEED output vector
261   Vec Eloc;
262   CeedVector eloc;
263   ierr = DMCreateLocalVector(dmEnergy, &Eloc); CHKERRQ(ierr);
264   ierr = VecGetSize(Eloc, &length); CHKERRQ(ierr);
265   ierr = VecDestroy(&Eloc); CHKERRQ(ierr);
266   CeedVectorCreate(user->ceed, length, &eloc);
267 
268   // Apply libCEED operator
269   CeedOperatorApply(opEnergy, user->Xceed, eloc, CEED_REQUEST_IMMEDIATE);
270 
271   // Restore PETSc vector
272   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
273   CHKERRQ(ierr);
274 
275   // Reduce max error
276   const CeedScalar *e;
277   CeedVectorGetArrayRead(eloc, CEED_MEM_HOST, &e);
278   (*energy) = 0;
279   for (CeedInt i=0; i<length; i++)
280     (*energy) += e[i];
281   CeedVectorRestoreArrayRead(eloc, &e);
282   CeedVectorDestroy(&eloc);
283 
284   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
285                        user->comm); CHKERRQ(ierr);
286 
287   PetscFunctionReturn(0);
288 };
289