Actual source code: petscda.h

  1: #pragma once

  3: #include <petscmat.h>

  5: /* MANSEC = ML */
  6: /* SUBMANSEC = PetscDA */

  8: /*S
  9:    PetscDA - Abstract PETSc object that manages data assimilation

 11:    Level: beginner

 13:    Notes:
 14:    This is new code, please independently verify all results you obtain using it.

 16:    Some planned work for `PetscDA` is available as GitLab Issue #1882

 18:    Currently we supply two ensemble-based assimilators: `PETSCDAETKF` and `PETSCDALETKF`

 20: .seealso: [](ch_da), `PetscDAType`, `PETSCDAETKF`, `PETSCDALETKF`, `PetscDASqrtType`, `PetscDACreate()`, `PetscDASetType()`,
 21:           `PetscDASetSizes()`, `PetscDAEnsembleSetSize()`, `PetscDAEnsembleAnalysis()`, `PetscDAEnsembleForecast()`,
 22:           `PetscDADestroy()`, `PetscDAView()`
 23: S*/
 24: typedef struct _p_PetscDA *PetscDA;

 26: /*E
 27:   PetscDASqrtType - Type of square root of matrices to use the data assimilation algorithms

 29:   Values:
 30: +  `PETSCDA_SQRT_CHOLESKY` - Use the Cholesky factorization
 31: -  `PETSCDA_SQRT_EIGEN`    - Use the eigenvalue decomposition

 33:   Option Database Key:
 34: . -petscda_ensemble_sqrt_type <cholesky, eigen> - select the square root type at run time

 36:   Level: intermediate

 38: .seealso: [](ch_da), `PetscDA`, `PetscDAEnsembleSetSqrtType()`, `PetscDAEnsembleGetSqrtType()`
 39: E*/
 40: typedef enum {
 41:   PETSCDA_SQRT_CHOLESKY = 0,
 42:   PETSCDA_SQRT_EIGEN    = 1
 43: } PetscDASqrtType;

 45: /*J
 46:   PetscDAType - String with the name of a PETSc data assimilation method

 48:   Level: beginner

 50: .seealso: [](ch_da), `PetscDA`, `PetscDASetType()`, `PETSCDAETKF`, `PETSCDALETKF`
 51: J*/
 52: typedef const char *PetscDAType;
 53: #define PETSCDAETKF  "etkf"
 54: #define PETSCDALETKF "letkf"

 56: PETSC_EXTERN PetscErrorCode PetscDAInitializePackage(void);
 57: PETSC_EXTERN PetscErrorCode PetscDAFinalizePackage(void);

 59: PETSC_EXTERN PetscErrorCode PetscDARegister(const char[], PetscErrorCode (*)(PetscDA));
 60: PETSC_EXTERN PetscErrorCode PetscDARegisterAll(void);

 62: PETSC_EXTERN PetscErrorCode PetscDACreate(MPI_Comm, PetscDA *);
 63: PETSC_EXTERN PetscErrorCode PetscDADestroy(PetscDA *);
 64: PETSC_EXTERN PetscErrorCode PetscDASetType(PetscDA, PetscDAType);
 65: PETSC_EXTERN PetscErrorCode PetscDAGetType(PetscDA, PetscDAType *);
 66: PETSC_EXTERN PetscErrorCode PetscDAView(PetscDA, PetscViewer);
 67: PETSC_EXTERN PetscErrorCode PetscDAViewFromOptions(PetscDA, PetscObject, const char[]);
 68: PETSC_EXTERN PetscErrorCode PetscDASetFromOptions(PetscDA);

 70: PETSC_EXTERN PetscErrorCode PetscDASetSizes(PetscDA, PetscInt, PetscInt);
 71: PETSC_EXTERN PetscErrorCode PetscDASetLocalSizes(PetscDA, PetscInt, PetscInt);
 72: PETSC_EXTERN PetscErrorCode PetscDAGetSizes(PetscDA, PetscInt *, PetscInt *);
 73: PETSC_EXTERN PetscErrorCode PetscDASetNDOF(PetscDA, PetscInt);
 74: PETSC_EXTERN PetscErrorCode PetscDAGetNDOF(PetscDA, PetscInt *);
 75: PETSC_EXTERN PetscErrorCode PetscDASetUp(PetscDA);

 77: PETSC_EXTERN PetscErrorCode PetscDASetObsErrorVariance(PetscDA, Vec);
 78: PETSC_EXTERN PetscErrorCode PetscDAGetObsErrorVariance(PetscDA, Vec *);

 80: PETSC_EXTERN PetscErrorCode PetscDASetOptionsPrefix(PetscDA, const char[]);
 81: PETSC_EXTERN PetscErrorCode PetscDAAppendOptionsPrefix(PetscDA, const char[]);
 82: PETSC_EXTERN PetscErrorCode PetscDAGetOptionsPrefix(PetscDA, const char *[]);

 84: PETSC_EXTERN PetscErrorCode PetscDAEnsembleSetSize(PetscDA, PetscInt);
 85: PETSC_EXTERN PetscErrorCode PetscDAEnsembleGetSize(PetscDA, PetscInt *);
 86: PETSC_EXTERN PetscErrorCode PetscDAEnsembleSetInflation(PetscDA, PetscReal);
 87: PETSC_EXTERN PetscErrorCode PetscDAEnsembleGetInflation(PetscDA, PetscReal *);

 89: PETSC_EXTERN PetscErrorCode PetscDAEnsembleGetMember(PetscDA, PetscInt, Vec *);
 90: PETSC_EXTERN PetscErrorCode PetscDAEnsembleRestoreMember(PetscDA, PetscInt, Vec *);
 91: PETSC_EXTERN PetscErrorCode PetscDAEnsembleSetMember(PetscDA, PetscInt, Vec);

 93: PETSC_EXTERN PetscErrorCode PetscDAEnsembleComputeMean(PetscDA, Vec);
 94: PETSC_EXTERN PetscErrorCode PetscDAEnsembleComputeAnomalies(PetscDA, Vec, Mat *);
 95: PETSC_EXTERN PetscErrorCode PetscDAEnsembleAnalysis(PetscDA, Vec, Mat);
 96: PETSC_EXTERN PetscErrorCode PetscDAEnsembleForecast(PetscDA, PetscErrorCode (*)(Vec, Vec, PetscCtx), PetscCtx);
 97: PETSC_EXTERN PetscErrorCode PetscDAEnsembleInitialize(PetscDA, Vec, PetscReal, PetscRandom);

 99: PETSC_EXTERN PetscErrorCode PetscDAEnsembleComputeNormalizedInnovationMatrix(Mat, Vec, Vec, PetscInt, PetscScalar, Mat);
100: PETSC_EXTERN PetscErrorCode PetscDAEnsembleSetSqrtType(PetscDA, PetscDASqrtType);
101: PETSC_EXTERN PetscErrorCode PetscDAEnsembleGetSqrtType(PetscDA, PetscDASqrtType *);
102: PETSC_EXTERN PetscErrorCode PetscDAEnsembleTFactor(PetscDA, Mat);
103: PETSC_EXTERN PetscErrorCode PetscDAEnsembleApplyTInverse(PetscDA, Vec, Vec);
104: PETSC_EXTERN PetscErrorCode PetscDAEnsembleApplySqrtTInverse(PetscDA, Mat, Mat);

106: PETSC_EXTERN PetscErrorCode PetscDALETKFSetLocalization(PetscDA, Mat, Mat);
107: PETSC_EXTERN PetscErrorCode PetscDALETKFSetObsPerVertex(PetscDA, PetscInt);
108: PETSC_EXTERN PetscErrorCode PetscDALETKFGetObsPerVertex(PetscDA, PetscInt *);

110: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
111: PETSC_EXTERN PetscErrorCode PetscDALETKFGetLocalizationMatrix(const PetscInt, const PetscInt, Vec[3], PetscReal[3], Mat, Mat *);
112: #endif