Actual source code: taotermutils.c
1: #include <petsc/private/taoimpl.h>
3: PETSC_INTERN PetscErrorCode VecIfNotCongruentGetSameLayoutVec(Vec a, Vec *b)
4: {
5: PetscFunctionBegin;
6: if (*b == NULL) {
7: PetscCall(VecDuplicate(a, b));
8: } else {
9: PetscLayout layout_a, layout_b;
10: PetscBool is_same;
12: PetscCall(VecGetLayout(a, &layout_a));
13: PetscCall(VecGetLayout(*b, &layout_b));
14: PetscCall(PetscLayoutCompare(layout_a, layout_b, &is_same));
15: if (is_same == PETSC_FALSE) {
16: PetscCall(VecDestroy(b));
17: PetscCall(VecDuplicate(a, b));
18: }
19: }
20: PetscFunctionReturn(PETSC_SUCCESS);
21: }
23: PETSC_INTERN PetscErrorCode TaoTermCreateHessianMatricesDefault_H_Internal(TaoTerm term, Mat *H, Mat *Hpre, PetscBool Hpre_is_H, MatType H_mattype)
24: {
25: Mat Htemp;
26: PetscBool is_mffd = PETSC_FALSE;
28: PetscFunctionBegin;
29: PetscCall(PetscStrcmp(H_mattype, MATMFFD, &is_mffd));
30: if (is_mffd) {
31: PetscCall(TaoTermCreateHessianMFFD(term, &Htemp));
32: } else {
33: PetscLayout sol_layout;
34: VecType sol_vec_type;
36: PetscCall(MatCreate(PetscObjectComm((PetscObject)term), &Htemp));
37: PetscCall(TaoTermGetSolutionLayout(term, &sol_layout));
38: PetscCall(MatSetLayouts(Htemp, sol_layout, sol_layout));
39: PetscCall(TaoTermGetSolutionVecType(term, &sol_vec_type));
40: if (H_mattype) PetscCall(MatSetType(Htemp, H_mattype));
41: else PetscCall(MatSetVecType(Htemp, sol_vec_type));
42: PetscCall(MatSetOption(Htemp, MAT_SYMMETRIC, PETSC_TRUE));
43: PetscCall(MatSetOption(Htemp, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
44: }
46: if (H) {
47: PetscCall(PetscObjectReference((PetscObject)Htemp));
48: *H = Htemp;
49: }
50: if (Hpre && Hpre_is_H) {
51: PetscCall(PetscObjectReference((PetscObject)Htemp));
52: *Hpre = Htemp;
53: }
54: PetscCall(MatDestroy(&Htemp));
55: PetscFunctionReturn(PETSC_SUCCESS);
56: }
58: PETSC_INTERN PetscErrorCode TaoTermCreateHessianMatricesDefault_Hpre_Internal(TaoTerm term, Mat *H, Mat *Hpre, PetscBool Hpre_is_H, MatType Hpre_mattype)
59: {
60: Mat Hpretemp;
61: PetscBool is_mffd = PETSC_FALSE;
63: PetscFunctionBegin;
64: PetscCall(PetscStrcmp(Hpre_mattype, MATMFFD, &is_mffd));
65: if (is_mffd) {
66: PetscCall(TaoTermCreateHessianMFFD(term, &Hpretemp));
67: } else {
68: PetscLayout sol_layout;
69: VecType sol_vec_type;
71: PetscCall(MatCreate(PetscObjectComm((PetscObject)term), &Hpretemp));
72: PetscCall(TaoTermGetSolutionLayout(term, &sol_layout));
73: PetscCall(MatSetLayouts(Hpretemp, sol_layout, sol_layout));
74: PetscCall(TaoTermGetSolutionVecType(term, &sol_vec_type));
75: if (Hpre_mattype) PetscCall(MatSetType(Hpretemp, Hpre_mattype));
76: else PetscCall(MatSetVecType(Hpretemp, sol_vec_type));
77: PetscCall(MatSetOption(Hpretemp, MAT_SYMMETRIC, PETSC_TRUE));
78: PetscCall(MatSetOption(Hpretemp, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
79: }
80: *Hpre = Hpretemp;
81: PetscFunctionReturn(PETSC_SUCCESS);
82: }