xref: /libCEED/examples/solids/src/matops.c (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*3d8e8822SJeremy L Thompson //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5*3d8e8822SJeremy L Thompson //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7*3d8e8822SJeremy L Thompson 
8ccaff030SJeremy L Thompson /// @file
9ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc
10ccaff030SJeremy L Thompson 
115754ecacSJeremy L Thompson #include "../include/matops.h"
125754ecacSJeremy L Thompson #include "../include/structs.h"
135754ecacSJeremy L Thompson #include "../include/utils.h"
14ccaff030SJeremy L Thompson 
15ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
16ccaff030SJeremy L Thompson // libCEED Operators for MatShell
17ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
18ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
19ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
20ccaff030SJeremy L Thompson   PetscErrorCode ierr;
21ccaff030SJeremy L Thompson   PetscScalar *x, *y;
22d1d35e2fSjeremylt   PetscMemType x_mem_type, y_mem_type;
23ccaff030SJeremy L Thompson 
24ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
25ccaff030SJeremy L Thompson 
26ccaff030SJeremy L Thompson   // Global-to-local
27d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
28d1d35e2fSjeremylt   ierr = VecZeroEntries(user->Y_loc); CHKERRQ(ierr);
29ccaff030SJeremy L Thompson 
30ccaff030SJeremy L Thompson   // Setup CEED vectors
31d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
32d1d35e2fSjeremylt                                    &x_mem_type); CHKERRQ(ierr);
33d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
34d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
35d1d35e2fSjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
36ccaff030SJeremy L Thompson 
37ccaff030SJeremy L Thompson   // Apply CEED operator
38d1d35e2fSjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, CEED_REQUEST_IMMEDIATE);
39ccaff030SJeremy L Thompson 
40ccaff030SJeremy L Thompson   // Restore PETSc vectors
41d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
42d1d35e2fSjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
43d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
44ccaff030SJeremy L Thompson   CHKERRQ(ierr);
45d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
46ccaff030SJeremy L Thompson 
47ccaff030SJeremy L Thompson   // Local-to-global
48ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
49d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, Y); CHKERRQ(ierr);
50ccaff030SJeremy L Thompson 
51ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
52ccaff030SJeremy L Thompson };
53ccaff030SJeremy L Thompson 
54ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
55ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
56ccaff030SJeremy L Thompson   PetscErrorCode ierr;
57ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
58ccaff030SJeremy L Thompson 
59ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
60ccaff030SJeremy L Thompson 
61ccaff030SJeremy L Thompson   // Use computed BCs
62d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
63d1d35e2fSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc,
64d1d35e2fSjeremylt                                     user->load_increment, NULL, NULL, NULL);
65ccaff030SJeremy L Thompson   CHKERRQ(ierr);
66ccaff030SJeremy L Thompson 
67ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
68ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
69ccaff030SJeremy L Thompson 
70fe394131Sjeremylt   // Neumann BCs
71d1d35e2fSjeremylt   if (user->neumann_bcs) {
72d1d35e2fSjeremylt     ierr = VecAXPY(Y, -user->load_increment, user->neumann_bcs); CHKERRQ(ierr);
73fe394131Sjeremylt   }
74fe394131Sjeremylt 
75ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
76ccaff030SJeremy L Thompson };
77ccaff030SJeremy L Thompson 
78ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
79ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
80ccaff030SJeremy L Thompson   PetscErrorCode ierr;
81ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
82ccaff030SJeremy L Thompson 
83ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
84ccaff030SJeremy L Thompson 
8562e9c006SJeremy L Thompson   // Zero boundary values
86d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
87ccaff030SJeremy L Thompson 
88ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
89ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
90ccaff030SJeremy L Thompson 
91ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
92ccaff030SJeremy L Thompson };
93ccaff030SJeremy L Thompson 
94ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
95ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
96ccaff030SJeremy L Thompson   PetscErrorCode ierr;
97ccaff030SJeremy L Thompson   UserMult user;
98ccaff030SJeremy L Thompson 
99ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
100ccaff030SJeremy L Thompson 
101ccaff030SJeremy L Thompson   // Zero boundary values
102ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
103d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
104ccaff030SJeremy L Thompson 
105ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
106ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
107ccaff030SJeremy L Thompson 
108ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
109ccaff030SJeremy L Thompson };
110ccaff030SJeremy L Thompson 
111ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
112ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
113ccaff030SJeremy L Thompson   PetscErrorCode ierr;
114ccaff030SJeremy L Thompson   UserMultProlongRestr user;
115ccaff030SJeremy L Thompson   PetscScalar *c, *f;
116d1d35e2fSjeremylt   PetscMemType c_mem_type, f_mem_type;
117ccaff030SJeremy L Thompson 
118ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
119ccaff030SJeremy L Thompson 
120ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
121ccaff030SJeremy L Thompson 
122ccaff030SJeremy L Thompson   // Global-to-local
123d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr);
124d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm_c, X, INSERT_VALUES, user->loc_vec_c);
125ccaff030SJeremy L Thompson   CHKERRQ(ierr);
126d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr);
127ccaff030SJeremy L Thompson 
128ccaff030SJeremy L Thompson   // Setup CEED vectors
129d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c,
130d1d35e2fSjeremylt                                    &c_mem_type); CHKERRQ(ierr);
131d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_f, &f, &f_mem_type); CHKERRQ(ierr);
132d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
133d1d35e2fSjeremylt                      c);
134d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
135d1d35e2fSjeremylt                      f);
136ccaff030SJeremy L Thompson 
137ccaff030SJeremy L Thompson   // Apply CEED operator
138d1d35e2fSjeremylt   CeedOperatorApply(user->op_prolong, user->ceed_vec_c, user->ceed_vec_f,
139ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
140ccaff030SJeremy L Thompson 
141ccaff030SJeremy L Thompson   // Restore PETSc vectors
142d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
143d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
144d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c);
145ccaff030SJeremy L Thompson   CHKERRQ(ierr);
146d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_f, &f); CHKERRQ(ierr);
147ccaff030SJeremy L Thompson 
148ccaff030SJeremy L Thompson   // Local-to-global
149ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
150d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm_f, user->loc_vec_f, ADD_VALUES, Y);
151ccaff030SJeremy L Thompson   CHKERRQ(ierr);
152ccaff030SJeremy L Thompson 
153ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
154ccaff030SJeremy L Thompson }
155ccaff030SJeremy L Thompson 
156ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
157ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
158ccaff030SJeremy L Thompson   PetscErrorCode ierr;
159ccaff030SJeremy L Thompson   UserMultProlongRestr user;
160ccaff030SJeremy L Thompson   PetscScalar *c, *f;
161d1d35e2fSjeremylt   PetscMemType c_mem_type, f_mem_type;
162ccaff030SJeremy L Thompson 
163ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
164ccaff030SJeremy L Thompson 
165ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
166ccaff030SJeremy L Thompson 
167ccaff030SJeremy L Thompson   // Global-to-local
168d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr);
169d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm_f, X, INSERT_VALUES, user->loc_vec_f);
170ccaff030SJeremy L Thompson   CHKERRQ(ierr);
171d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr);
172ccaff030SJeremy L Thompson 
173ccaff030SJeremy L Thompson   // Setup CEED vectors
174d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f,
175d1d35e2fSjeremylt                                    &f_mem_type); CHKERRQ(ierr);
176d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_c, &c, &c_mem_type); CHKERRQ(ierr);
177d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
178d1d35e2fSjeremylt                      f);
179d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
180d1d35e2fSjeremylt                      c);
181ccaff030SJeremy L Thompson 
182ccaff030SJeremy L Thompson   // Apply CEED operator
183d1d35e2fSjeremylt   CeedOperatorApply(user->op_restrict, user->ceed_vec_f, user->ceed_vec_c,
184ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
185ccaff030SJeremy L Thompson 
186ccaff030SJeremy L Thompson   // Restore PETSc vectors
187d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
188d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
189d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f);
190ccaff030SJeremy L Thompson   CHKERRQ(ierr);
191d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_c, &c); CHKERRQ(ierr);
192ccaff030SJeremy L Thompson 
193ccaff030SJeremy L Thompson   // Local-to-global
194ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
195d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm_c, user->loc_vec_c, ADD_VALUES, Y);
196ccaff030SJeremy L Thompson   CHKERRQ(ierr);
197ccaff030SJeremy L Thompson 
198ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
199ccaff030SJeremy L Thompson };
2002d93065eSjeremylt 
201ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
202ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
203ccaff030SJeremy L Thompson   PetscErrorCode ierr;
204ccaff030SJeremy L Thompson   UserMult user;
205ccaff030SJeremy L Thompson 
206ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
207ccaff030SJeremy L Thompson 
208ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
209ccaff030SJeremy L Thompson 
210f7b4142eSJeremy L Thompson   // -- Set physics context
211d1d35e2fSjeremylt   if (user->ctx_phys_smoother)
212d1d35e2fSjeremylt     CeedQFunctionSetContext(user->qf, user->ctx_phys_smoother);
213f7b4142eSJeremy L Thompson 
2142bba3ffaSJeremy L Thompson   // Compute Diagonal via libCEED
2152bba3ffaSJeremy L Thompson   PetscScalar *x;
216d1d35e2fSjeremylt   PetscMemType x_mem_type;
2172bba3ffaSJeremy L Thompson 
2182bba3ffaSJeremy L Thompson   // -- Place PETSc vector in libCEED vector
219d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->X_loc, &x, &x_mem_type); CHKERRQ(ierr);
220d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2212bba3ffaSJeremy L Thompson 
2222bba3ffaSJeremy L Thompson   // -- Compute Diagonal
223d1d35e2fSjeremylt   CeedOperatorLinearAssembleDiagonal(user->op, user->x_ceed,
224ccaff030SJeremy L Thompson                                      CEED_REQUEST_IMMEDIATE);
225ccaff030SJeremy L Thompson 
226f7b4142eSJeremy L Thompson   // -- Reset physics context
227d1d35e2fSjeremylt   if (user->ctx_phys_smoother)
228d1d35e2fSjeremylt     CeedQFunctionSetContext(user->qf, user->ctx_phys);
229f7b4142eSJeremy L Thompson 
230ccaff030SJeremy L Thompson   // -- Local-to-Global
231d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
232d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->X_loc, &x); CHKERRQ(ierr);
233ccaff030SJeremy L Thompson   ierr = VecZeroEntries(D); CHKERRQ(ierr);
234d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm, user->X_loc, ADD_VALUES, D); CHKERRQ(ierr);
235ccaff030SJeremy L Thompson 
2362bba3ffaSJeremy L Thompson   // Cleanup
237d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
238ccaff030SJeremy L Thompson 
239ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
240ccaff030SJeremy L Thompson };
2412d93065eSjeremylt 
2422d93065eSjeremylt // This function calculates the strain energy in the final solution
243a3c02c40SJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
244d1d35e2fSjeremylt                                    CeedOperator op_energy, Vec X,
245a3c02c40SJeremy L Thompson                                    PetscReal *energy) {
2462d93065eSjeremylt   PetscErrorCode ierr;
2472d93065eSjeremylt   PetscScalar *x;
248d1d35e2fSjeremylt   PetscMemType x_mem_type;
2492d93065eSjeremylt   CeedInt length;
2502d93065eSjeremylt 
2512d93065eSjeremylt   PetscFunctionBeginUser;
2522d93065eSjeremylt 
2532d93065eSjeremylt   // Global-to-local
254d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
255d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
256d1d35e2fSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc,
257d1d35e2fSjeremylt                                     user->load_increment, NULL, NULL, NULL);
2585c25879aSJeremy L Thompson   CHKERRQ(ierr);
2592d93065eSjeremylt 
260a3c02c40SJeremy L Thompson   // Setup libCEED input vector
261d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
262d1d35e2fSjeremylt                                    &x_mem_type); CHKERRQ(ierr);
263d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2642d93065eSjeremylt 
265a3c02c40SJeremy L Thompson   // Setup libCEED output vector
266d1d35e2fSjeremylt   Vec E_loc;
267d1d35e2fSjeremylt   CeedVector e_loc;
268d1d35e2fSjeremylt   ierr = DMCreateLocalVector(dmEnergy, &E_loc); CHKERRQ(ierr);
269d1d35e2fSjeremylt   ierr = VecGetSize(E_loc, &length); CHKERRQ(ierr);
270d1d35e2fSjeremylt   ierr = VecDestroy(&E_loc); CHKERRQ(ierr);
271d1d35e2fSjeremylt   CeedVectorCreate(user->ceed, length, &e_loc);
272a3c02c40SJeremy L Thompson 
2732d93065eSjeremylt   // Apply libCEED operator
274d1d35e2fSjeremylt   CeedOperatorApply(op_energy, user->x_ceed, e_loc, CEED_REQUEST_IMMEDIATE);
2752d93065eSjeremylt 
2762d93065eSjeremylt   // Restore PETSc vector
277d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
278d1d35e2fSjeremylt   ierr = VecRestoreArrayRead(user->X_loc, (const PetscScalar **)&x);
2792d93065eSjeremylt   CHKERRQ(ierr);
2802d93065eSjeremylt 
2812d93065eSjeremylt   // Reduce max error
2822d93065eSjeremylt   const CeedScalar *e;
283d1d35e2fSjeremylt   CeedVectorGetArrayRead(e_loc, CEED_MEM_HOST, &e);
2842d93065eSjeremylt   (*energy) = 0;
2852d93065eSjeremylt   for (CeedInt i=0; i<length; i++)
2862d93065eSjeremylt     (*energy) += e[i];
287d1d35e2fSjeremylt   CeedVectorRestoreArrayRead(e_loc, &e);
288d1d35e2fSjeremylt   CeedVectorDestroy(&e_loc);
2892d93065eSjeremylt 
2902d93065eSjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
2912d93065eSjeremylt                        user->comm); CHKERRQ(ierr);
2922d93065eSjeremylt 
2932d93065eSjeremylt   PetscFunctionReturn(0);
2942d93065eSjeremylt };
295