Developer Notes¶
Style Guide¶
Please check your code for style issues by running
make format
which uses clang-format
.
All style fixes provided by make format
should be committed.
In addition to those automatically enforced style rules, Ratel tends to follow the following code style conventions:
Variable names:
snake_case
Struct members:
snake_case
Function names:
PascalCase
Type names:
PascalCase
Constant names:
CAPS_SNAKE_CASE
Also, documentation files should have one sentence per line to help make git diffs clearer and less disruptive.
Clang-tidy¶
Please check your code for common issues by running
make tidy
which uses the clang-tidy
utility included in recent releases of Clang.
This tool is much slower than actual compilation (make -j
parallelism helps).
All issues reported by make tidy
should be fixed.
Header Files¶
Header inclusion for source files should follow the principal of ‘include what you use’ rather than relying upon transitive #include
to define all symbols.
Every symbol that is used in the source file foo.c
should be defined in foo.c
, foo.h
, or in a header file #include
d in one of these two locations.
Please check your code by running the tool include-what-you-use
to see recommendations for changes to your source.
Most issues reported by include-what-you-use
should be fixed; however this rule is flexible to account for differences in header file organization in external libraries.
Header files should be listed in alphabetical order, with installed headers preceding local headers.
#include <ceed.h>
#include <petsc.h>
#include <stdbool.h>
#include <string.h>
#include "include/ratel.h"
restrict
Semantics¶
QFunction arguments can be assumed to have restrict
semantics.
That is, each input and output array must reside in distinct memory without overlap.
Continuous Integration (CI) and Testing¶
Whenever a merge request is created, the Ratel test suite is run against a variety of hardware and software configurations.
If one of the tests in the CI pipeline fails, you should look at the error in the CI job log.
The log will provide the command which failed and relevant details (such as diff
output) to help debug the issue.
You can also see the CGNS or CSV output of a failed run, if applicable, by selecting “Download” on the Job Artifacts (to the right of the log).
Artifacts from failed runs will be contained in the test_failure_artifacts
folder inside the artifacts archive.
Adding New Boundary Conditions¶
Each boundary condition (e.g. pressure, slip, traction, flux, etc.) gets its own folder include/ratel/boundary/[bc-name]
containing:
params.h
: Necessary data structures, must be JiT-safeqf.h
: QFunctions for setup (optional) and action(s) of the boundary condition, (obviously) must be JiT-safe[bc-name].h
: Header for functions specific to the boundary condition, which will in the future be standardized into a uniform interface. This header should include./params.h
, but NOT include./qf.h
.
An example for [bc-name].h
for non-Dirichlet BCs is include/ratel/boundary/traction/traction.h
:
/// @file
/// Ratel Traction Boundary Conditions Functions
#pragma once
#include <ratel-impl.h>
#include "params.h"
RATEL_INTERN PetscErrorCode RatelBoundaryTractionDataFromOptions(Ratel ratel, PetscInt i, RatelBCTractionParams *params_traction);
RATEL_INTERN PetscErrorCode RatelBoundaryTractionSetupEnergySuboperators(Ratel ratel, DM dm_energy, CeedOperator op_external_energy);
RATEL_INTERN PetscErrorCode RatelBoundaryTractionSetupSuboperators(Ratel ratel, DM dm, CeedOperator op_residual);
For Dirichlet BCs, see include/ratel/boundary/clamp/clamp.h
:
/// @file
/// Ratel Clamp Boundary Conditions Functions
#pragma once
#include <ratel-impl.h>
#include "params.h"
RATEL_INTERN PetscErrorCode RatelBoundaryClampDataFromOptions(Ratel ratel, PetscInt i, PetscInt j, RatelBCClampParams *params_clamp);
RATEL_INTERN PetscErrorCode RatelBoundaryClampSetupDirichletSuboperators(Ratel ratel, DM dm, PetscBool incremental, CeedOperator op_dirichlet);
The functions that should exist in the header are a subset of the following, along with any needed helper functions:
RatelBoundary[bc-name]DataFromOptions
, reads command-line options into the relevant data structure(s) defined inparams.h
, may also take integers forface_id
and/orfield_index
RatelBoundary[bc-name]SetupDirichletSuboperators
, for Dirichlet boundariesRatelBoundary[bc-name]SetupSuboperators
, for adding suboperators to any ofop_residual_u
,op_residual_ut
,op_residual_utt
,op_jacobian
RatelBoundary[bc-name]SetupMaterialSuboperators
, like*SetupSuboperators
, but for BCs with explicit material dependence (right now only Nitsche contact). The first argument should beRatelMaterial
instead ofRatel
.RatelBoundary[bc-name]SetupEnergySuboperators
, for adding suboperators to the external energy operator, if applicable
These will eventually be changed to a unified interface under a RatelBoundary
data structure.
The implementations for these headers go in the corresponding src/boundary/[bc-name].c
file.
The existing header ratel-boundary.h
now contains utility functions for setting up boundary conditions.
Most critical of these are:
RatelBoundarySetupSuboperators
: all BCs which define aRatelBoundary[bc-name]SetupSuboperators
should have those called here. If a new suboperator may be added toop_jacobian
, special care must be taken to set up the multigrid level (see below).RatelBoundarySetupMaterialSuboperators
: all BCs which define aRatelBoundary[bc-name]SetupMaterialSuboperators
should be called here. See below for details onop_jacobian
suboperators.RatelBoundarySetupDirichletSuboperators
: all BCs which define aRatelBoundary[bc-name]SetupDirichletSuboperators
should be called hereRatelBoundarySetupEnergySuboperators
: all BCs which define aRatelBoundary[bc-name]SetupEnergySuboperators
should be called here
BCs with Jacobian suboperators¶
There is now a simple process for setting up multigrid levels for BCs with Jacobian suboperators.
Now, both penalty contact and pressure BCs use RatelMaterialSetBoundaryJacobianMultigridInfo
and the new RatelBoundarySetupMultigridLevel_Standard
and RatelBoundarySetupMultigridLevel_CellToFace
helpers to ensure that MG levels are set up for all Jacobian suboperators.
As an implementation detail, the Jacobian indices for these boundaries are considered to live on the last material, but this choice is arbitrary and not relevant to the process of adding new boundary conditions.
For non-material dependent BCs, the call to RatelBoundary[bc-name]SetupSuboperators
should be wrapped like:
const RatelMaterial last_material = ratel->materials[ratel->num_materials - 1];
CeedInt old_num_jacobian_sub_operators, new_num_jacobian_sub_operators;
RatelCallCeed(ratel, CeedCompositeOperatorGetNumSub(op_jacobian, &old_num_jacobian_sub_operators));
PetscCall(RatelBoundary[bc-name]SetupSuboperators(ratel, dm, op_residual_u, op_jacobian));
RatelCallCeed(ratel, CeedCompositeOperatorGetNumSub(op_jacobian, &new_num_jacobian_sub_operators));
// -- Set Jacobian info on last material
PetscCall(RatelMaterialSetBoundaryJacobianMultigridInfo(last_material, op_jacobian, old_num_jacobian_sub_operators, new_num_jacobian_sub_operators,
&RatelBoundarySetupMultigridLevel_Standard));
For material-dependent BCs, the call to RatelBoundary[bc-name]SetupMaterialSuboperators
should be wrapped like:
RatelCallCeed(ratel, CeedCompositeOperatorGetNumSub(op_jacobian, &old_num_jacobian_sub_operators));
PetscCall(RatelBoundary[bc-name]SetupMaterialSuboperators(material, u_dot, op_residual_u, op_jacobian));
RatelCallCeed(ratel, CeedCompositeOperatorGetNumSub(op_jacobian, &new_num_jacobian_sub_operators));
// Set Jacobian info for material
PetscCall(RatelMaterialSetBoundaryJacobianMultigridInfo(material, op_jacobian, old_num_jacobian_sub_operators, new_num_jacobian_sub_operators,
&RatelBoundarySetupMultigridLevel_CellToFace));
Adding New Initial Conditions¶
Ratel supports non-zero initial condition for linear poromechanincs MMS and as a constant value for other cases. If user needs to add non-zero initial condition as function of coordinate, a function with the same structure as RatelSetupInitialConditionMMS
needs to be added with its own Qfucntion expressing the initial condition of each active field.