PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
DirichletData.hh
Go to the documentation of this file.
1 /* Copyright (C) 2020, 2022 PISM Authors
2  *
3  * This file is part of PISM.
4  *
5  * PISM is free software; you can redistribute it and/or modify it under the
6  * terms of the GNU General Public License as published by the Free Software
7  * Foundation; either version 3 of the License, or (at your option) any later
8  * version.
9  *
10  * PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11  * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13  * details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with PISM; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 #ifndef PISM_DIRICHLETDATA_H
20 #define PISM_DIRICHLETDATA_H
21 
22 #include "pism/util/fem/FEM.hh"
23 #include "pism/util/petscwrappers/Mat.hh" // Mat
24 #include "pism/util/Vector2d.hh"
25 
26 namespace pism {
27 
28 namespace array {
29 class Array;
30 class Scalar;
31 class Vector;
32 }
33 
34 
35 namespace fem {
36 
37 class Element2;
38 
39 //* Parts shared by scalar and 2D vector Dirichlet data classes.
41 public:
42  void constrain(Element2 &element);
43  operator bool() {
44  return m_indices != NULL;
45  }
46 protected:
47  DirichletData();
49 
50  void init(const array::Scalar *indices, const array::Array *values, double weight = 1.0);
51  void finish(const array::Array *values);
52 
55  double m_weight;
56 };
57 
59 public:
60  DirichletData_Scalar(const array::Scalar *indices, const array::Scalar *values,
61  double weight = 1.0);
63 
64  void enforce(const Element2 &element, double* x_e);
65  void enforce_homogeneous(const Element2 &element, double* x_e);
66  void fix_residual(double const *const *const x_global, double **r_global);
67  void fix_residual_homogeneous(double **r_global);
68  void fix_jacobian(Mat J);
69 protected:
71 };
72 
74 public:
75  DirichletData_Vector(const array::Scalar *indices, const array::Vector *values,
76  double weight);
78 
79  void enforce(const Element2 &element, Vector2d* x_e);
80  void enforce_homogeneous(const Element2 &element, Vector2d* x_e);
81  void fix_residual(Vector2d const *const *const x_global, Vector2d **r_global);
83  void fix_jacobian(Mat J);
84 protected:
86 };
87 
88 } // end of namespace fem
89 } // end of namespace pism
90 
91 #endif /* PISM_DIRICHLETDATA_H */
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition: Vector2d.hh:29
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition: Array.hh:208
void enforce_homogeneous(const Element2 &element, double *x_e)
void fix_residual(double const *const *const x_global, double **r_global)
void fix_residual_homogeneous(double **r_global)
const array::Scalar * m_values
DirichletData_Scalar(const array::Scalar *indices, const array::Scalar *values, double weight=1.0)
void enforce(const Element2 &element, double *x_e)
void fix_residual(Vector2d const *const *const x_global, Vector2d **r_global)
DirichletData_Vector(const array::Scalar *indices, const array::Vector *values, double weight)
void enforce(const Element2 &element, Vector2d *x_e)
void fix_residual_homogeneous(Vector2d **r)
void enforce_homogeneous(const Element2 &element, Vector2d *x_e)
const array::Vector * m_values
double m_indices_e[q1::n_chi]
void finish(const array::Array *values)
void init(const array::Scalar *indices, const array::Array *values, double weight=1.0)
const array::Scalar * m_indices
void constrain(Element2 &element)
Constrain element, i.e. ensure that quadratures do not contribute to Dirichlet nodes by marking corre...
const int n_chi
Definition: FEM.hh:191
static int weight(int M_ij, int M_n, double h_ij, double h_n)
const int J[]
Definition: ssafd_code.cc:34