xref: /libCEED/examples/solids/src/matops.c (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33d8e8822SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53d8e8822SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d8e8822SJeremy 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"
122b730f8bSJeremy L Thompson 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdmplex.h>
1549aac155SJeremy L Thompson 
165754ecacSJeremy L Thompson #include "../include/structs.h"
175754ecacSJeremy L Thompson #include "../include/utils.h"
18ccaff030SJeremy L Thompson 
19ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
20ccaff030SJeremy L Thompson // libCEED Operators for MatShell
21ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
22ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
ApplyLocalCeedOp(Vec X,Vec Y,UserMult user)23ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
24ccaff030SJeremy L Thompson   PetscScalar *x, *y;
25d1d35e2fSjeremylt   PetscMemType x_mem_type, y_mem_type;
26ccaff030SJeremy L Thompson 
27ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
28ccaff030SJeremy L Thompson 
29ccaff030SJeremy L Thompson   // Global-to-local
302b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc));
312b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->Y_loc));
32ccaff030SJeremy L Thompson 
33ccaff030SJeremy L Thompson   // Setup CEED vectors
342b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, &x_mem_type));
352b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type));
36d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
37d1d35e2fSjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
38ccaff030SJeremy L Thompson 
39ccaff030SJeremy L Thompson   // Apply CEED operator
40d1d35e2fSjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, CEED_REQUEST_IMMEDIATE);
41ccaff030SJeremy L Thompson 
42ccaff030SJeremy L Thompson   // Restore PETSc vectors
43d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
44d1d35e2fSjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
452b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x));
462b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(user->Y_loc, &y));
47ccaff030SJeremy L Thompson 
48ccaff030SJeremy L Thompson   // Local-to-global
492b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
502b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, Y));
51ccaff030SJeremy L Thompson 
52ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53ccaff030SJeremy L Thompson };
54ccaff030SJeremy L Thompson 
55ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
FormResidual_Ceed(SNES snes,Vec X,Vec Y,void * ctx)56ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
57ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
58ccaff030SJeremy L Thompson 
59ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
60ccaff030SJeremy L Thompson 
61ccaff030SJeremy L Thompson   // Use computed BCs
622b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->X_loc));
632b730f8bSJeremy L Thompson   PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc, user->load_increment, NULL, NULL, NULL));
64ccaff030SJeremy L Thompson 
65ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
662b730f8bSJeremy L Thompson   PetscCall(ApplyLocalCeedOp(X, Y, user));
67ccaff030SJeremy L Thompson 
68fe394131Sjeremylt   // Neumann BCs
69d1d35e2fSjeremylt   if (user->neumann_bcs) {
702b730f8bSJeremy L Thompson     PetscCall(VecAXPY(Y, -user->load_increment, user->neumann_bcs));
71fe394131Sjeremylt   }
72fe394131Sjeremylt 
73ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
74ccaff030SJeremy L Thompson };
75ccaff030SJeremy L Thompson 
76ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
ApplyJacobianCoarse_Ceed(SNES snes,Vec X,Vec Y,void * ctx)77ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
78ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
79ccaff030SJeremy L Thompson 
80ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
81ccaff030SJeremy L Thompson 
8262e9c006SJeremy L Thompson   // Zero boundary values
832b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->X_loc));
84ccaff030SJeremy L Thompson 
85ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
862b730f8bSJeremy L Thompson   PetscCall(ApplyLocalCeedOp(X, Y, user));
87ccaff030SJeremy L Thompson 
88ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
89ccaff030SJeremy L Thompson };
90ccaff030SJeremy L Thompson 
91ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
ApplyJacobian_Ceed(Mat A,Vec X,Vec Y)92ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
93ccaff030SJeremy L Thompson   UserMult user;
94ccaff030SJeremy L Thompson 
95ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
96ccaff030SJeremy L Thompson 
97ccaff030SJeremy L Thompson   // Zero boundary values
982b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &user));
992b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->X_loc));
100ccaff030SJeremy L Thompson 
101ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
1022b730f8bSJeremy L Thompson   PetscCall(ApplyLocalCeedOp(X, Y, user));
103ccaff030SJeremy L Thompson 
104ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
105ccaff030SJeremy L Thompson };
106ccaff030SJeremy L Thompson 
107ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
Prolong_Ceed(Mat A,Vec X,Vec Y)108ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
109ccaff030SJeremy L Thompson   UserMultProlongRestr user;
110ccaff030SJeremy L Thompson   PetscScalar         *c, *f;
111d1d35e2fSjeremylt   PetscMemType         c_mem_type, f_mem_type;
112ccaff030SJeremy L Thompson 
113ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
114ccaff030SJeremy L Thompson 
1152b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &user));
116ccaff030SJeremy L Thompson 
117ccaff030SJeremy L Thompson   // Global-to-local
1182b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->loc_vec_c));
1192b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm_c, X, INSERT_VALUES, user->loc_vec_c));
1202b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->loc_vec_f));
121ccaff030SJeremy L Thompson 
122ccaff030SJeremy L Thompson   // Setup CEED vectors
1232b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c, &c_mem_type));
1242b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(user->loc_vec_f, &f, &f_mem_type));
1252b730f8bSJeremy L Thompson   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER, c);
1262b730f8bSJeremy L Thompson   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER, f);
127ccaff030SJeremy L Thompson 
128ccaff030SJeremy L Thompson   // Apply CEED operator
1292b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_prolong, user->ceed_vec_c, user->ceed_vec_f, CEED_REQUEST_IMMEDIATE);
130ccaff030SJeremy L Thompson 
131ccaff030SJeremy L Thompson   // Restore PETSc vectors
132d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
133d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
1342b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c));
1352b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(user->loc_vec_f, &f));
136ccaff030SJeremy L Thompson 
137ccaff030SJeremy L Thompson   // Local-to-global
1382b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
1392b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm_f, user->loc_vec_f, ADD_VALUES, Y));
140ccaff030SJeremy L Thompson 
141ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
142ccaff030SJeremy L Thompson }
143ccaff030SJeremy L Thompson 
144ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
Restrict_Ceed(Mat A,Vec X,Vec Y)145ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
146ccaff030SJeremy L Thompson   UserMultProlongRestr user;
147ccaff030SJeremy L Thompson   PetscScalar         *c, *f;
148d1d35e2fSjeremylt   PetscMemType         c_mem_type, f_mem_type;
149ccaff030SJeremy L Thompson 
150ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
151ccaff030SJeremy L Thompson 
1522b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &user));
153ccaff030SJeremy L Thompson 
154ccaff030SJeremy L Thompson   // Global-to-local
1552b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->loc_vec_f));
1562b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm_f, X, INSERT_VALUES, user->loc_vec_f));
1572b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->loc_vec_c));
158ccaff030SJeremy L Thompson 
159ccaff030SJeremy L Thompson   // Setup CEED vectors
1602b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f, &f_mem_type));
1612b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(user->loc_vec_c, &c, &c_mem_type));
1622b730f8bSJeremy L Thompson   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER, f);
1632b730f8bSJeremy L Thompson   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER, c);
164ccaff030SJeremy L Thompson 
165ccaff030SJeremy L Thompson   // Apply CEED operator
1662b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_restrict, user->ceed_vec_f, user->ceed_vec_c, CEED_REQUEST_IMMEDIATE);
167ccaff030SJeremy L Thompson 
168ccaff030SJeremy L Thompson   // Restore PETSc vectors
169d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
170d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
1712b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f));
1722b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(user->loc_vec_c, &c));
173ccaff030SJeremy L Thompson 
174ccaff030SJeremy L Thompson   // Local-to-global
1752b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Y));
1762b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm_c, user->loc_vec_c, ADD_VALUES, Y));
177ccaff030SJeremy L Thompson 
178ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
179ccaff030SJeremy L Thompson };
1802d93065eSjeremylt 
181ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
GetDiag_Ceed(Mat A,Vec D)182ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
183ccaff030SJeremy L Thompson   UserMult user;
184ccaff030SJeremy L Thompson 
185ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
186ccaff030SJeremy L Thompson 
1872b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(A, &user));
188ccaff030SJeremy L Thompson 
189f7b4142eSJeremy L Thompson   // -- Set physics context
1902b730f8bSJeremy L Thompson   if (user->ctx_phys_smoother) CeedQFunctionSetContext(user->qf, user->ctx_phys_smoother);
191f7b4142eSJeremy L Thompson 
1922bba3ffaSJeremy L Thompson   // Compute Diagonal via libCEED
1932bba3ffaSJeremy L Thompson   PetscScalar *x;
194d1d35e2fSjeremylt   PetscMemType x_mem_type;
1952bba3ffaSJeremy L Thompson 
1962bba3ffaSJeremy L Thompson   // -- Place PETSc vector in libCEED vector
1972b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(user->X_loc, &x, &x_mem_type));
198d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
1992bba3ffaSJeremy L Thompson 
2002bba3ffaSJeremy L Thompson   // -- Compute Diagonal
2012b730f8bSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op, user->x_ceed, CEED_REQUEST_IMMEDIATE);
202ccaff030SJeremy L Thompson 
203f7b4142eSJeremy L Thompson   // -- Reset physics context
2042b730f8bSJeremy L Thompson   if (user->ctx_phys_smoother) CeedQFunctionSetContext(user->qf, user->ctx_phys);
205f7b4142eSJeremy L Thompson 
206ccaff030SJeremy L Thompson   // -- Local-to-Global
207d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
2082b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(user->X_loc, &x));
2092b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(D));
2102b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, user->X_loc, ADD_VALUES, D));
211ccaff030SJeremy L Thompson 
2122bba3ffaSJeremy L Thompson   // Cleanup
2132b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->X_loc));
214ccaff030SJeremy L Thompson 
215ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
216ccaff030SJeremy L Thompson };
2172d93065eSjeremylt 
2182d93065eSjeremylt // This function calculates the strain energy in the final solution
ComputeStrainEnergy(DM dmEnergy,UserMult user,CeedOperator op_energy,Vec X,PetscReal * energy)2192b730f8bSJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user, CeedOperator op_energy, Vec X, PetscReal *energy) {
2202d93065eSjeremylt   PetscScalar *x;
221d1d35e2fSjeremylt   PetscMemType x_mem_type;
2224d00b080SJeremy L Thompson   PetscInt     length;
2232d93065eSjeremylt 
2242d93065eSjeremylt   PetscFunctionBeginUser;
2252d93065eSjeremylt 
2262d93065eSjeremylt   // Global-to-local
2272b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(user->X_loc));
2282b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc));
2292b730f8bSJeremy L Thompson   PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc, user->load_increment, NULL, NULL, NULL));
2302d93065eSjeremylt 
231a3c02c40SJeremy L Thompson   // Setup libCEED input vector
2322b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, &x_mem_type));
233d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2342d93065eSjeremylt 
235a3c02c40SJeremy L Thompson   // Setup libCEED output vector
236d1d35e2fSjeremylt   Vec        E_loc;
237d1d35e2fSjeremylt   CeedVector e_loc;
2384d00b080SJeremy L Thompson 
2392b730f8bSJeremy L Thompson   PetscCall(DMCreateLocalVector(dmEnergy, &E_loc));
2402b730f8bSJeremy L Thompson   PetscCall(VecGetSize(E_loc, &length));
2412b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&E_loc));
242d1d35e2fSjeremylt   CeedVectorCreate(user->ceed, length, &e_loc);
243a3c02c40SJeremy L Thompson 
2442d93065eSjeremylt   // Apply libCEED operator
245d1d35e2fSjeremylt   CeedOperatorApply(op_energy, user->x_ceed, e_loc, CEED_REQUEST_IMMEDIATE);
2462d93065eSjeremylt 
2472d93065eSjeremylt   // Restore PETSc vector
248d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
2492b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayRead(user->X_loc, (const PetscScalar **)&x));
2502d93065eSjeremylt 
2512d93065eSjeremylt   // Reduce max error
2522d93065eSjeremylt   const CeedScalar *e;
253d1d35e2fSjeremylt   CeedVectorGetArrayRead(e_loc, CEED_MEM_HOST, &e);
2542d93065eSjeremylt   (*energy) = 0;
2552b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < length; i++) (*energy) += e[i];
256d1d35e2fSjeremylt   CeedVectorRestoreArrayRead(e_loc, &e);
257d1d35e2fSjeremylt   CeedVectorDestroy(&e_loc);
2582d93065eSjeremylt 
2592b730f8bSJeremy L Thompson   PetscCall(MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM, user->comm));
2602d93065eSjeremylt 
261ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2622d93065eSjeremylt };
263